openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
electrics.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 
26 #include "electrics.h"
27 #include "petsc_utils.h"
28 #include "timers.h"
29 #include "stimulate.h"
30 #include "electric_integrators.h"
31 
32 #include "SF_init.h" // for SF::init_xxx()
33 
34 #ifdef WITH_CALIPER
35 #include "caliper/cali.h"
36 #else
37 #include "caliper_hooks.h"
38 #endif
39 
40 
41 namespace opencarp {
42 
44 {
46  double t1, t2;
47  get_time(t1);
48 
49  set_dir(OUTPUT);
50 
51  // open logger
52  logger = f_open("electrics.log", param_globals::experiment != 4 ? "w" : "r");
53 
54  // setup mappings between extra and intra grids, algebraic and nodal,
55  // and between PETSc and canonical orderings
56  setup_mappings();
57 
58  // the ionic physics is currently triggered from inside the Electrics to have tighter
59  // control over it
60  ion.logger = logger;
61  ion.initialize();
62 
63  // set up Intracellular tissue
65  region_mask(intra_elec_msh, mtype[intra_grid].regions, mtype[intra_grid].regionIDs, true, "gregion_i");
66 
67  if (param_globals::bidomain || param_globals::extracell_monodomain_stim) {
68  // set up Extracellular tissue
70  region_mask(extra_elec_msh, mtype[extra_grid].regions, mtype[extra_grid].regionIDs, true, "gregion_e");
71  }
72 
73  // add electrics timer for time stepping, add to time stepper tool (TS)
74  double global_time = user_globals::tm_manager->time;
75  timer_idx = user_globals::tm_manager->add_eq_timer(global_time, param_globals::tend, 0,
76  param_globals::dt, 0, "elec::ref_dt", "TS");
77 
78  // electrics stimuli setup
79  setup_stimuli();
80 
81  // set up the linear equation systems. this needs to happen after the stimuli have been
82  // set up, since we need boundary condition info
83  setup_solvers();
84 
85  // the next setup steps require the solvers to be set up, since they use the matrices
86  // generated by those
87 
88  // balance electrodes, we may need the extracellular mass matrix
89  balance_electrodes();
90  // total current scaling
92  // initialize the LATs detector
94 
95  // initialize phie recovery data
96  if(strlen(param_globals::phie_rec_ptf) > 0)
98 
99  // prepare the electrics output. we skip it if we do post-processing
100  if(param_globals::experiment != EXP_POSTPROCESS)
101  setup_output();
102 
103  if (param_globals::prepacing_bcl > 0)
104  prepace();
105  this->initialize_time += timing(t2, t1);
106 }
107 
109 {
110  MaterialType *m = mtype+g;
111 
112  // initialize random conductivity fluctuation structure with PrM values
113  m->regions.resize(param_globals::num_gregions);
114 
115  const char* grid_name = g == Electrics::intra_grid ? "intracellular" : "extracellular";
116  log_msg(logger, 0, 0, "Setting up %s tissue properties for %d regions ..", grid_name,
117  param_globals::num_gregions);
118 
119  char buf[64];
120  RegionSpecs* reg = m->regions.data();
121 
122  for (size_t i=0; i<m->regions.size(); i++, reg++) {
123  if(!strcmp(param_globals::gregion[i].name, "")) {
124  snprintf(buf, sizeof buf, ", gregion_%d", int(i));
125  param_globals::gregion[i].name = dupstr(buf);
126  }
127 
128  reg->regname = strdup(param_globals::gregion[i].name);
129  reg->regID = i;
130  reg->nsubregs = param_globals::gregion[i].num_IDs;
131  if(!reg->nsubregs)
132  reg->subregtags = NULL;
133  else {
134  reg->subregtags = new int[reg->nsubregs];
135  for (int j=0;j<reg->nsubregs;j++) {
136  reg->subregtags[j] = param_globals::gregion[i].ID[j];
137  if(reg->subregtags[j]==-1)
138  log_msg(NULL,3,ECHO, "Warning: not all %u IDs provided for gregion[%u]!\n", reg->nsubregs, i);
139  }
140  }
141 
142  // describe material in given region
143  elecMaterial *emat = new elecMaterial();
144  emat->material_type = ElecMat;
145 
146  emat->InVal[0] = param_globals::gregion[i].g_il;
147  emat->InVal[1] = param_globals::gregion[i].g_it;
148  emat->InVal[2] = param_globals::gregion[i].g_in;
149 
150  emat->ExVal[0] = param_globals::gregion[i].g_el;
151  emat->ExVal[1] = param_globals::gregion[i].g_et;
152  emat->ExVal[2] = param_globals::gregion[i].g_en;
153 
154  emat->BathVal[0] = param_globals::gregion[i].g_bath;
155  emat->BathVal[1] = param_globals::gregion[i].g_bath;
156  emat->BathVal[2] = param_globals::gregion[i].g_bath;
157 
158  // convert units from S/m -> mS/um
159  for (int j=0; j<3; j++) {
160  emat->InVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
161  emat->ExVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
162  emat->BathVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
163  }
164  reg->material = emat;
165  }
166 
167  {
169  const char* file = g == Electrics::intra_grid ? param_globals::gi_scale_vec : param_globals::ge_scale_vec;
170  if (strlen(file))
171  read_el_scale_vec(file, mt, m->el_scale, m->el_scale_dpn);
172  }
173 }
174 
175 void Electrics::setup_mappings()
176 {
178  bool intra_exits = mesh_is_registered(intra_elec_msh), extra_exists = mesh_is_registered(extra_elec_msh);
179  assert(intra_exits);
180  const int dpn = 1;
181 
182  // It may be that another physic (e.g. ionic models) has already computed the intracellular mappings,
183  // thus we first test their existence
184  if(get_scattering(intra_elec_msh, ALG_TO_NODAL, dpn) == NULL) {
185  log_msg(logger, 0, 0, "%s: Setting up intracellular algebraic-to-nodal scattering.", __func__);
187  }
189  log_msg(logger, 0, 0, "%s: Setting up intracellular PETSc to canonical permutation.", __func__);
191  }
192 
193  // extracellular mappings
194  if(extra_exists) {
195  log_msg(logger, 0, 0, "%s: Setting up extracellular algebraic-to-nodal scattering.", __func__);
197  log_msg(logger, 0, 0, "%s: Setting up extracellular PETSc to canonical permutation.", __func__);
199  log_msg(logger, 0, 0, "%s: Setting up intra-to-extra scattering.", __func__);
201  }
202 
203  bool check_i2e = false;
204  if(check_i2e && extra_exists) {
205  sf_mesh & intra_mesh = get_mesh(intra_elec_msh);
206  sf_mesh & extra_mesh = get_mesh(extra_elec_msh);
207  int rank = get_rank();
208 
210 
211  const SF::vector<mesh_int_t> & intra_alg_nod = intra_mesh.pl.algebraic_nodes();
212  const SF::vector<mesh_int_t> & extra_alg_nod = extra_mesh.pl.algebraic_nodes();
213  const SF::vector<mesh_int_t> & extra_petsc_nbr = extra_mesh.get_numbering(SF::NBR_PETSC);
214  const SF::vector<mesh_int_t> & intra_ref_nbr = intra_mesh.get_numbering(SF::NBR_REF);
215  const SF::vector<mesh_int_t> & extra_ref_nbr = extra_mesh.get_numbering(SF::NBR_REF);
216 
217  // TODO(init) : delete these three at the end of this section?
218  sf_vec *intra_testvec; SF::init_vector(&intra_testvec, intra_mesh, 1, sf_vec::algebraic);
219  sf_vec *extra_testvec; SF::init_vector(&extra_testvec, extra_mesh, 1, sf_vec::algebraic);
220  sf_vec *i2e_testvec; SF::init_vector(&i2e_testvec, extra_mesh, 1, sf_vec::algebraic);
221 
222  SF_real* id = intra_testvec->ptr();
223  for(size_t i=0; i<intra_alg_nod.size(); i++) {
224  int lpidx = local_nodal_to_local_petsc(intra_mesh, rank, intra_alg_nod[i]);
225  id[lpidx] = intra_ref_nbr[intra_alg_nod[i]];
226  }
227  intra_testvec->release_ptr(id);
228 
229  SF_real* ed = extra_testvec->ptr();
230  for(size_t i=0; i<extra_alg_nod.size(); i++) {
231  int lpidx = local_nodal_to_local_petsc(extra_mesh, rank, extra_alg_nod[i]);
232  ed[lpidx] = extra_ref_nbr[extra_alg_nod[i]];
233  }
234  extra_testvec->release_ptr(ed);
235 
236  i2e_testvec->set(-1.0);
237  i2e.forward(*intra_testvec, *i2e_testvec);
238 
239  int err = 0;
240  for(size_t i=0; i<extra_alg_nod.size(); i++) {
241  auto id = i2e_testvec->get(i);
242  auto ed = extra_testvec->get(i);
243  if(id > -1 && id != ed)
244  err++;
245  }
246 
247  if(get_global(err, MPI_SUM))
248  log_msg(0,5,0, "Electrics mapping test failed!");
249  else
250  log_msg(0,5,0, "Electrics mapping test succeeded!");
251  }
252 }
253 
255 {
257  double t1, t2;
258  get_time(t1);
259 
260  // if requested, we checkpoint the current state
261  checkpointing();
262 
263  // activation checking
264  const double time = user_globals::tm_manager->time,
265  time_step = user_globals::tm_manager->time_step;
266  lat.check_acts(time);
267  lat.check_quiescence(time, time_step);
268 
269  // I believe that we need to treat the stimuli in two ways:
270  // - Extracellular potential stimuli (this includes ground) affect the
271  // elliptic solver in a more delicate way, as such, there is a dbc_manager
272  // to take care of that.
273  // - Extracellular and Intracellular current stimuli are applied to the rhs vectors
274  // and can be managed by the stimulate() code directly.
275  stimulate_extracellular();
276 
277 if(param_globals::bidomain == BIDOMAIN)
279 
280  clamp_Vm();
281 
282  // compute ionics update
283  ion.compute_step();
284 
285  stimulate_intracellular();
286 
287  // store Vm before parabolic step, the full Ic we compute in the output step
288  if(param_globals::dump_data & DUMP_IC)
290 
291  // solver parabolic system
293 
294  clamp_Vm();
295 
296  if(user_globals::tm_manager->trigger(iotm_console)) {
297  // output lin solver stats
299  if(param_globals::bidomain == BIDOMAIN)
301  }
302  this->compute_time += timing(t2, t1);
303 
304  // since the traces have their own timing, we check for trace dumps in the compute step loop
307 
308 }
309 
311 {
313  double t1, t2;
314  get_time(t1);
315 
316  const double time = user_globals::tm_manager->time,
317  time_step = user_globals::tm_manager->time_step;
318 
319  // for pseudo-bidomain we compute extracellular potential only for output
320  if(param_globals::bidomain == PSEUDO_BIDM) {
323  ellip_solver.stats.log_stats(time, false);
324  }
325 
326  if(param_globals::dump_data & DUMP_IVOL)
328 
329  if(param_globals::bidomain && (param_globals::dump_data & DUMP_IACT)) {
331  }
332 
333  if(param_globals::dump_data & DUMP_IC) {
334  PetscReal *Ic = parab_solver.Ic->ptr(), *Vmv = parab_solver.Vmv->ptr();
335 
336  if (time <= time_step) {
337  for(PetscInt i=0; i < parab_solver.Ic->lsize(); i++)
338  Ic[i] = (Vmv[i] - Vmv[i]) / (-time_step); // first timestep: finite difference = 0 (no previous Vmv available)
339  }
340  else{
341  for(PetscInt i=0; i < parab_solver.Ic->lsize(); i++)
342  Ic[i] = (Ic[i] - Vmv[i]) / (-time_step);
343  }
344 
345 
347  }
348 
349  // recover phie
350  if(phie_rcv.pts.size()) {
352  }
353 
356 
357  double curtime = timing(t2, t1);
358  this->output_time += curtime;
359 
360  IO_stats.calls++;
361  IO_stats.tot_time += curtime;
362 
364  IO_stats.log_stats(time, false);
365 }
366 
371 {
373  // output LAT data
375 
376  // close logger
377  f_close(logger);
378 
379  // close output files
381 
382  // destroy ionics
383  ion.destroy();
384 }
385 
386 void balance_electrode(elliptic_solver & ellip, SF::vector<stimulus> & stimuli, int balance_from, int balance_to)
387 {
388  log_msg( NULL, 0, 0, "Balancing stimulus %d with %d %s-wise.",balance_from, balance_to,
389  is_current(stimuli[balance_from].phys.type) ? "current" : "voltage" );
390 
391  stimulus & from = stimuli[balance_from];
392  stimulus & to = stimuli[balance_to];
393 
394  to.pulse = from.pulse;
395  to.ptcl = from.ptcl;
396  to.phys = from.phys;
397  to.pulse.strength *= -1.0;
398 
399  if (from.phys.type == I_ex)
400  {
401  // if from is total current, skip volume based adjustment of strength
402  // otherwise, calling constant_total_stimulus_current() will undo the balanced scaling of to.pulse.strength
403  // constant_total_stimulus_current() will do the scaling based on the volume
404  if (!from.phys.total_current) {
405  sf_mat& mass = *ellip.mass_e;
406  SF_real vol0 = get_volume_from_nodes(mass, from.electrode.vertices);
408 
409  to.pulse.strength *= fabs(vol0 / vol1);
410  }
411  }
412 }
413 
414 void Electrics::balance_electrodes()
415 {
416  for(int i=0; i<param_globals::num_stim; i++) {
417  if(param_globals::stim[i].crct.balance != -1) {
418  int from = param_globals::stim[i].crct.balance;
419  int to = i;
420 
421  balance_electrode(this->ellip_solver, stimuli, from, to);
422  }
423  }
424 }
425 
426 void Electrics::setup_stimuli()
427 {
428  // initialize basic stim info data (used units, supported types, etc)
429  init_stim_info();
430 
431  stimuli.resize(param_globals::num_stim);
432  for(int i=0; i<param_globals::num_stim; i++)
433  {
434  // construct new stimulus
435  stimulus & s = stimuli[i];
436 
438  s.translate(i);
439 
440  s.setup(i);
441 
442  if (s.electrode.dump_vtx)
443  s.dump_vtx_file(i);
444 
445  if(param_globals::stim[i].pulse.dumpTrace && get_rank() == 0) {
446  set_dir(OUTPUT);
447  s.pulse.wave.write_trace(s.name+".trc");
448  }
449  }
450 }
451 
452 void apply_stim_to_vector(const stimulus & s, sf_vec & vec, bool add)
453 {
454  double val; s.value(val);
455  const SF::vector<mesh_int_t> & idx = s.electrode.vertices;
456  const int rank = get_rank();
457  SF::vector<mesh_int_t> local_idx = idx;
458  for (size_t i = 0; i < idx.size(); i++) {
459  local_idx[i] = local_nodal_to_local_petsc(*vec.mesh, rank, idx[i]);
460  }
461  vec.set(local_idx, val, add, true);
462 }
463 
464 void Electrics::stimulate_intracellular()
465 {
466  parabolic_solver & ps = parab_solver;
467 
468  // iterate over stimuli
469  for(stimulus & s : stimuli) {
470  if(s.is_active()) {
471  // for active stimuli, deal with the stimuli-type specific stimulus application
472  switch(s.phys.type)
473  {
474  case I_tm: {
475  if(param_globals::operator_splitting) {
476  apply_stim_to_vector(s, *ps.Vmv, true);
477  }
478  else {
479  SF_real Cm = 1.0;
480  timer_manager & tm = *user_globals::tm_manager;
481  SF_real sc = tm.time_step / Cm;
482 
483  ps.Irhs->set(0.0);
484  apply_stim_to_vector(s, *ps.Irhs, true);
485 
486  *ps.tmp_i1 = *ps.IIon;
487  *ps.tmp_i1 -= *ps.Irhs;
488  *ps.tmp_i1 *= sc; // tmp_i1 = sc * (IIon - Irhs)
489 
490  // add ionic, transmembrane and intracellular currents to rhs
491  if(param_globals::parab_solve != parabolic_solver::EXPLICIT)
492  ps.mass_i->mult(*ps.tmp_i1, *ps.Irhs);
493  else
494  *ps.Irhs = *ps.tmp_i1;
495  }
496  break;
497  }
498 
499  case Illum: {
500  sf_vec* illum_vec = ion.miif->gdata[limpet::illum];
501 
502  if(illum_vec == NULL) {
503  log_msg(0,5,0, "Cannot apply illumination stim: global vector not present!");
504  EXIT(EXIT_FAILURE);
505  } else {
506  apply_stim_to_vector(s, *illum_vec, false);
507  }
508 
509  break;
510  }
511 
512  default: break;
513  }
514  }
515  }
516 }
517 
518 void Electrics::clamp_Vm() {
519  for(stimulus & s : stimuli) {
520  if(s.phys.type == Vm_clmp && s.is_active())
522  }
523 }
524 
525 void Electrics::stimulate_extracellular()
526 {
527  if(param_globals::bidomain) {
528  // we check if the DBC layout changed, if so we recompute the matrix and the dbc_manager
529  bool dbcs_have_updated = ellip_solver.dbc != nullptr && ellip_solver.dbc->dbc_update();
531 
532  if(dbcs_have_updated && time_not_final)
534 
535  ellip_solver.phiesrc->set(0.0);
536 
537  for(const stimulus & s : stimuli) {
538  if(s.is_active() && s.phys.type == I_ex)
540  }
541  }
542 }
543 
545  SF::vector<mesh_int_t> & inp_idx,
547 {
548  int mpi_rank = get_rank(), mpi_size = get_size();
549  const SF::vector<mesh_int_t> & layout = mesh.pl.algebraic_layout();
550 
551  SF::vector<mesh_int_t> sndbuff;
552 
553  size_t buffsize = 0;
554  idx.resize(0);
555 
556  for(int pid=0; pid < mpi_size; pid++) {
557  if(mpi_rank == pid) {
558  sndbuff = inp_idx;
559  buffsize = sndbuff.size();
560  }
561 
562  MPI_Bcast(&buffsize, sizeof(size_t), MPI_BYTE, pid, PETSC_COMM_WORLD);
563  sndbuff.resize(buffsize);
564  MPI_Bcast(sndbuff.data(), buffsize*sizeof(mesh_int_t), MPI_BYTE, pid, PETSC_COMM_WORLD);
565 
566  mesh_int_t start = layout[mpi_rank], stop = layout[mpi_rank+1];
567 
568  for(mesh_int_t i : sndbuff) {
569  if(i >= start && i < stop)
570  idx.push_back(i - start);
571  }
572  }
573 
574  binary_sort(idx); unique_resize(idx);
575 }
576 
578  SF::vector<mesh_int_t> & inp_idx,
580 {
581  int mpi_rank = get_rank(), mpi_size = get_size();
582  const SF::vector<mesh_int_t> & alg_nod = mesh.pl.algebraic_nodes();
584 
586  for(mesh_int_t ii : alg_nod)
587  amap[nbr[ii]] = ii;
588 
589  SF::vector<mesh_int_t> sndbuff;
590  size_t buffsize = 0;
591  idx.resize(0);
592 
593  for(int pid=0; pid < mpi_size; pid++) {
594  if(mpi_rank == pid) {
595  sndbuff = inp_idx;
596  buffsize = sndbuff.size();
597  }
598 
599  MPI_Bcast(&buffsize, sizeof(size_t), MPI_BYTE, pid, PETSC_COMM_WORLD);
600  sndbuff.resize(buffsize);
601  MPI_Bcast(sndbuff.data(), buffsize*sizeof(mesh_int_t), MPI_BYTE, pid, PETSC_COMM_WORLD);
602 
603  for(mesh_int_t i : sndbuff) {
604  if(amap.count(i))
605  idx.push_back(amap[i]);
606  }
607  }
608 
609  binary_sort(idx); unique_resize(idx);
610 }
611 
612 void setup_dataout(const int dataout, std::string dataout_vtx, mesh_t grid,
613  SF::vector<mesh_int_t>* & restr, bool async)
614 {
615  sf_mesh & mesh = get_mesh(grid);
616 
617  switch(dataout) {
618 
619  case DATAOUT_SURF: {
620  sf_mesh surfmesh;
621  compute_surface_mesh(mesh, SF::NBR_SUBMESH, surfmesh);
622 
623  SF::vector<mesh_int_t> idxbuff(surfmesh.con);
624  binary_sort(idxbuff); unique_resize(idxbuff);
625 
626  restr = new SF::vector<mesh_int_t>();
627 
628  // for sync output, we need restr to hold the local indices in the petsc vectors
629  // that have been permuted to canonical numbering. For async, we need the
630  // non-overlapping decomposition of indices in NBR_SUBMESH numbering. The petsc indices will be
631  // computed at a later stage. The only reason we need to call compute_restr_idx_async,
632  // is that surface nodes in NBR_SUBMESH, may
633  // reside on partitions where they are not part of the algebraic nodes. thus we need to
634  // recommunicate to make sure the data layout is correct. We do not have this problem for
635  // DATAOUT_VTX.
636  if(!async)
637  compute_restr_idx(mesh, idxbuff, *restr);
638  else
639  compute_restr_idx_async(mesh, idxbuff, *restr);
640 
641  break;
642  }
643 
644  case DATAOUT_VTX: {
645  SF::vector<mesh_int_t> idxbuff;
646 
647  update_cwd();
648 
649  set_dir(INPUT);
650  read_indices(idxbuff, dataout_vtx, mesh, SF::NBR_REF, true, PETSC_COMM_WORLD);
651  set_dir(CURDIR);
652 
653  restr = new SF::vector<mesh_int_t>();
654 
655  if(!async) {
657  for(mesh_int_t & i : idxbuff) i = nbr[i];
658 
659  compute_restr_idx(mesh, idxbuff, *restr);
660  } else {
661  *restr = idxbuff;
662  }
663 
664  break;
665  }
666 
667  case DATAOUT_NONE:
668  case DATAOUT_VOL:
669  default: break;
670  }
671 }
672 
673 void Electrics::setup_output()
674 {
675  int rank = get_rank();
676  SF::vector<mesh_int_t>* restr_i = NULL;
677  SF::vector<mesh_int_t>* restr_e = NULL;
678  set_dir(OUTPUT);
679 
680  setup_dataout(param_globals::dataout_i, param_globals::dataout_i_vtx, intra_elec_msh,
681  restr_i, param_globals::num_io_nodes > 0);
682 
683  if(param_globals::dataout_i)
684  output_manager.register_output(parab_solver.Vmv, intra_elec_msh, 1, param_globals::vofile, "mV", restr_i);
685 
686  if(param_globals::bidomain) {
687  setup_dataout(param_globals::dataout_e, param_globals::dataout_e_vtx, extra_elec_msh,
688  restr_e, param_globals::num_io_nodes > 0);
689 
690  if(param_globals::dataout_i)
691  output_manager.register_output(ellip_solver.phie_i, intra_elec_msh, 1, param_globals::phieifile, "mV", restr_i);
692  if(param_globals::dataout_e)
693  output_manager.register_output(ellip_solver.phie, extra_elec_msh, 1, param_globals::phiefile, "mV", restr_e);
694  }
695 
696  if(param_globals::dump_data & DUMP_IC) {
697  output_manager.register_output(parab_solver.Ic, intra_elec_msh, 1, "Ic.igb", "uA/cm^2", restr_i);
698  output_manager.register_output(parab_solver.IIon, intra_elec_msh, 1, "Iion.igb","uA/cm^2", restr_i);
699  }
700 
701  if(param_globals::dump_data & DUMP_IVOL)
702  output_manager.register_output(parab_solver.Ivol, intra_elec_msh, 1, "Ivol.igb", "uA", restr_i);
703  if(param_globals::dump_data & DUMP_IACT)
704  output_manager.register_output(parab_solver.Iact, intra_elec_msh, 1, "Iact.igb", "uA", restr_i);
705 
706  if(phie_rcv.pts.size())
707  output_manager.register_output_sync(phie_rcv.phie_rec, phie_recv_msh, 1, param_globals::phie_recovery_file, "mV");
708 
710 
711  if(param_globals::num_trace) {
712  sf_mesh & imesh = get_mesh(intra_elec_msh);
713  open_trace(ion.miif, param_globals::num_trace, param_globals::trace_node, NULL, &imesh);
714  }
715 
716  // initialize generic logger for IO timings per time_dt
717  IO_stats.init_logger("IO_stats.dat");
718 }
719 
720 void Electrics::dump_matrices()
721 {
722  std::string bsname = param_globals::dump_basename;
723  std::string fn;
724 
725  set_dir(OUTPUT);
726 
727  // dump monodomain matrices
728  if ( param_globals::parab_solve==1 ) {
729  // using Crank-Nicolson
730  fn = bsname + "_Ki_CN.bin";
731  parab_solver.lhs_parab->write(fn.c_str());
732  }
733  fn = bsname + "_Ki.bin";
734  parab_solver.rhs_parab->write(fn.c_str());
735 
736  fn = bsname + "_Mi.bin";
737  parab_solver.mass_i->write(fn.c_str());
738 
739  if ( param_globals::bidomain ) {
740  fn = bsname + "_Kie.bin";
741  ellip_solver.phie_mat->write(fn.c_str());
742 
743  fn = bsname + "_Me.bin";
744  ellip_solver.mass_e->write(fn.c_str());
745  }
746 }
747 
748 
751 double Electrics::timer_val(const int timer_id)
752 {
753  // determine
754  int sidx = stimidx_from_timeridx(stimuli, timer_id);
755  double val = 0.0;
756  if(sidx != -1) {
757  stimuli[sidx].value(val);
758  }
759  else
760  val = std::nan("NaN");
761 
762  return val;
763 }
764 
767 std::string Electrics::timer_unit(const int timer_id)
768 {
769  int sidx = stimidx_from_timeridx(stimuli, timer_id);
770  std::string s_unit;
771 
772  if(sidx != -1)
773  // found a timer-linked stimulus
774  s_unit = stimuli[sidx].pulse.wave.f_unit;
775 
776  return s_unit;
777 }
778 
781 int stimidx_from_timeridx(const SF::vector<stimulus> & stimuli, const int timer_id)
782 {
783  // the only electrical quantities linked to a timer are stimuli
784  // thus we search for timer links only among stimuli for now
785 
786  // iterate over stimuli
787  for(size_t i = 0; i<stimuli.size(); i++)
788  {
789  const stimulus & s = stimuli[i];
790 
791  if(s.ptcl.timer_id == timer_id)
792  return s.idx;
793  }
794 
795  // invalid timer index not linked to any stimulus
796  return -1;
797 }
798 
809 void get_kappa(sf_vec & kappa, IMPregion *ir, limpet::MULTI_IF & miif, double k)
810 {
811  double* reg_kappa = new double[miif.N_IIF];
812 
813  for(int i=0; i<miif.N_IIF; i++)
814  reg_kappa[i] = k * miif.IIF[i]->cgeom().SVratio * ir[i].volFrac;
815 
816  double *kd = kappa.ptr();
817 
818  for(int i = 0; i < miif.numNode; i++)
819  kd[i] = reg_kappa[(int) miif.IIFmask[i]];
820 
821  kappa.release_ptr(kd);
822  delete [] reg_kappa;
823 }
824 
825 
834 {
835  for(size_t i=0; i < m.regions.size(); i++) {
836  elecMaterial *emat = static_cast<elecMaterial*>(m.regions[i].material);
837  emat->g = type;
838  }
839 }
840 
841 void Electrics::setup_solvers()
842 {
843  set_dir(OUTPUT);
844  parab_solver.init();
846 
847  if (param_globals::bidomain) {
848  ellip_solver.init();
850  }
851 
852  if(param_globals::dump2MatLab)
853  dump_matrices();
854 }
855 
857 {
858  for(const stimulus & s : stimuli) {
859  if(is_dbc(s.phys.type))
860  return true;
861  }
862  return false;
863 }
864 
865 const char* get_tsav_ext(double time)
866 {
867  int min_idx = -1;
868  double min_diff = 1e100;
869 
870  for(int i=0; i<param_globals::num_tsav; i++)
871  {
872  double diff = fabs(param_globals::tsav[i] - time);
873  if(min_diff > diff) {
874  min_diff = diff;
875  min_idx = i;
876  }
877  }
878 
879  if(min_idx == -1)
880  min_idx = 0;
881 
882  return param_globals::tsav_ext[min_idx];
883 }
884 
885 void Electrics::checkpointing()
886 {
887  const timer_manager & tm = *user_globals::tm_manager;
888 
889  // regular user selected state save
890  if (tm.trigger(iotm_chkpt_list)) {
891  char save_fnm[1024];
892  const char* tsav_ext = get_tsav_ext(tm.time);
893 
894  snprintf(save_fnm, sizeof save_fnm, "%s.%s.roe", param_globals::write_statef, tsav_ext);
895 
896  ion.miif->dump_state(save_fnm, tm.time, intra_elec_msh, false, GIT_COMMIT_COUNT);
897  }
898 
899  // checkpointing based on interval
900  if (tm.trigger(iotm_chkpt_intv)) {
901  char save_fnm[1024];
902  snprintf(save_fnm, sizeof save_fnm, "checkpoint.%.1f.roe", tm.time);
903  ion.miif->dump_state(save_fnm, tm.time, intra_elec_msh, false, GIT_COMMIT_COUNT);
904  }
905 }
906 
908 {
910  double t0, t1, dur;
911  get_time(t0);
912  stats.init_logger("ell_stats.dat");
913 
914  // here we can differentiate the solvers
916  sf_mesh & extra_mesh = get_mesh(extra_elec_msh);
917  sf_vec::ltype alg_type = sf_vec::algebraic;
918  const int dpn = 1;
919 
920  SF::init_vector(&phie, extra_mesh, dpn, alg_type);
921  SF::init_vector(&phiesrc, extra_mesh, dpn, alg_type);
922  SF::init_vector(&currtmp, extra_mesh, dpn, alg_type);
923 
925  sf_mesh & intra_mesh = get_mesh(intra_elec_msh);
926  SF::init_vector(&phie_i, intra_mesh, dpn, alg_type);
927  }
928 
929  int max_row_entries = max_nodal_edgecount(extra_mesh);
930 
933 
934  // alloc stiffness matrix
935  phie_mat->init(extra_mesh, dpn, dpn, max_row_entries);
936  // alloc mass matrix
937  mass_e ->init(extra_mesh, dpn, dpn, param_globals::mass_lumping ? 1 : max_row_entries);
938  dur = timing(t1, t0);
939 }
940 
942  SF::vector<stimulus> & stimuli,
943  FILE_SPEC logger)
944 {
946  double t0, t1, dur;
947  get_time(t0);
948  rebuild_stiffness(mtype, stimuli, logger);
949  rebuild_mass(logger);
950  dur = timing(t1, t0);
951 }
952 
954  SF::vector<stimulus> & stimuli,
955  FILE_SPEC logger)
956 {
958  double t0, t1, dur;
959  int log_flag = param_globals::output_level > 1 ? ECHO : 0;
960 
961  MaterialType & mt = mtype[Electrics::extra_grid];
962  const bool have_dbc = have_dbc_stims(stimuli);
963 
964  cond_t condType = sum_cond;
965  set_cond_type(mt, condType);
966 
967  // get mesh reference
968  sf_mesh & mesh = get_mesh(extra_elec_msh);
969 
970  get_time(t0);
971 
972  // fill the system
973  elec_stiffness_integrator stfn_integ(mt);
974 
975  phie_mat->zero();
976  SF::assemble_matrix(*phie_mat, mesh, stfn_integ);
977  phie_mat->scale(-1.0);
978 
979  dur = timing(t1,t0);
980  log_msg(logger,0,log_flag, "Computed ellipitc stiffness matrix in %.3f seconds.", dur);
981 
982  // set boundary conditions
983  if(have_dbc) {
984  log_msg(logger,0,log_flag, "Elliptic lhs matrix enforcing Dirichlet boundaries.");
985  get_time(t0);
986 
987  if(dbc == nullptr)
988  dbc = new dbc_manager(*phie_mat, stimuli);
989  else
990  dbc->recompute_dbcs();
991 
992  dbc->enforce_dbc_lhs();
993 
994  dur = timing(t1,t0);
995  log_msg(logger,0,log_flag, "Elliptic lhs matrix Dirichlet enforcing done in %.3f seconds.", dur);
996  }
997  else {
998  log_msg(logger,1,ECHO, "Elliptic lhs matrix is singular!");
999  // we are dealing with a singular system
1000  phie_mat_has_nullspace = true;
1001  }
1002 
1003  // solver has not been initialized yet
1004  set_dir(INPUT);
1005  get_time(t0);
1006 
1007  setup_linear_solver(logger);
1008 
1009  dur = timing(t1,t0);
1010  log_msg(logger,0,log_flag, "Initializing elliptic solver in %.5f seconds.", dur);
1011  set_dir(OUTPUT);
1012 }
1013 
1015 {
1017  int log_flag = param_globals::output_level > 1 ? ECHO : 0;
1018  double t0, t1, dur;
1019  mass_integrator mass_integ;
1020 
1021  // get mesh reference
1022  sf_mesh & mesh = get_mesh(extra_elec_msh);
1023  get_time(t0);
1024  mass_e->zero();
1025 
1026  if(param_globals::mass_lumping) {
1027  SF::assemble_lumped_matrix(*mass_e, mesh, mass_integ);
1028  } else {
1029  SF::assemble_matrix(*mass_e, mesh, mass_integ);
1030  }
1031 
1032  dur = timing(t1,t0);
1033  log_msg(logger,0,log_flag, "Computed elliptic mass matrix in %.3f seconds.", dur);
1034 }
1035 
1036 void elliptic_solver::setup_linear_solver(FILE_SPEC logger)
1037 {
1039 
1040  tol = param_globals::cg_tol_ellip;
1041  max_it = param_globals::cg_maxit_ellip;
1042 
1043  std::string default_opts;
1044  std::string solver_file;
1045  solver_file = param_globals::ellip_options_file;
1046  if (param_globals::flavor == std::string("ginkgo")) {
1047  default_opts = std::string(
1048  R"({
1049  "type": "solver::Cg",
1050  "criteria": [
1051  {
1052  "type": "Iteration",
1053  "max_iters": 100
1054  },
1055  {
1056  "type": "ResidualNorm",
1057  "reduction_factor": 1e-4
1058  }
1059  ],
1060  "preconditioner": {
1061  "type": "solver::Multigrid",
1062  "mg_level": [
1063  {
1064  "type": "multigrid::Pgm",
1065  "deterministic": true
1066  }
1067  ],
1068  "criteria": [
1069  {
1070  "type": "Iteration",
1071  "max_iters": 1
1072  }
1073  ],
1074  "coarsest_solver": {
1075  "type": "preconditioner::Schwarz",
1076  "local_solver": {
1077  "type": "preconditioner::Ilu"
1078  }
1079  },
1080  "max_levels": 10,
1081  "min_coarse_rows": 8,
1082  "default_initial_guess": "zero"
1083  }
1084 })");
1085  } else if (param_globals::flavor == std::string("petsc")) {
1086  default_opts = std::string("-ksp_type cg -pc_type hypre -pc_hypre_type boomeramg -pc_hypre_boomeramg_max_iter 1 -pc_hypre_boomeramg_strong_threshold 0.0 -options_left");
1087  }
1088 
1089  lin_solver->setup_solver(*phie_mat, tol, max_it, param_globals::cg_norm_ellip,
1090  "elliptic PDE", phie_mat_has_nullspace,
1091  logger, solver_file.c_str(), default_opts.c_str());
1092 }
1093 
1094 void elliptic_solver::solve(sf_mat & Ki, sf_vec & Vmv, sf_vec & tmp_i)
1097  double t0,t1;
1099  // assembly of rhs for FE
1100  if (phiesrc->mag() > 0.0) {
1101  mass_e->mult(*phiesrc, *currtmp);
1103  }
1104 
1105  Ki.mult(Vmv, tmp_i);
1106 
1107  bool add = true;
1108  i2e->forward(tmp_i, *phiesrc, add);
1109 
1110  if(dbc != nullptr)
1112 
1113  get_time(t0);
1114  (*lin_solver)(*phie, *phiesrc);
1115 
1116  // treat solver statistics
1117  auto dur = timing(t1, t0);
1118  lin_solver->time += dur;
1119  stats.slvtime += dur;
1121  if(lin_solver->reason < 0) {
1122  log_msg(0, 5, 0,"%s solver diverged. Reason: %s.", lin_solver->name.c_str(),
1123  petsc_get_converged_reason_str(lin_solver->reason));
1124  EXIT(1);
1125  }
1127  add = false;
1128  i2e->backward(*phie, *phie_i, add);
1129 }
1130 
1132 {
1134  double t0,t1;
1135 
1136  if(dbc != nullptr)
1138 
1139  get_time(t0);
1140  (*lin_solver)(*phie, *phiesrc);
1141 
1142  // treat solver statistics
1143  auto dur = timing(t1, t0);
1144  lin_solver->time += dur;
1145  stats.slvtime += dur;
1147 
1148  if(lin_solver->reason < 0) {
1149  log_msg(0, 5, 0,"%s solver diverged. Reason: %s.", lin_solver->name.c_str(),
1150  petsc_get_converged_reason_str(lin_solver->reason));
1151  EXIT(1);
1152  }
1153 
1154  // phie_i is only set up when we have an IntraMesh registered
1155  if(is_init(phie_i)) {
1156  bool add = false;
1158  i2e->backward(*phie, *phie_i, add);
1159  }
1160 }
1161 
1163 {
1165  double t0, t1, dur;
1166  get_time(t0);
1167  stats.init_logger("par_stats.dat");
1168 
1169  // here we can differentiate the solvers
1171 
1172  sf_vec* vm_ptr = get_data(vm_vec);
1173  sf_vec* iion_ptr = get_data(iion_vec);
1174 
1175  if(!(vm_ptr != NULL && iion_ptr != NULL)) {
1176  log_msg(0,5,0, "%s error: global Vm and Iion vectors not properly set up! Ionics seem invalid! Aborting!",
1177  __func__);
1178  EXIT(1);
1179  }
1180 
1181  SF::init_vector(&Vmv);
1183  Vmv-> shallow_copy(*vm_ptr);
1184  IIon->shallow_copy(*iion_ptr);
1185 
1186  if(param_globals::dump_data & DUMP_IC) SF::init_vector(&Ic , Vmv);
1187  if(param_globals::dump_data & DUMP_IVOL) SF::init_vector(&Ivol, Vmv);
1188  if(param_globals::dump_data & DUMP_IACT) SF::init_vector(&Iact, Vmv);
1189 
1192 
1193  sf_mesh & intra_mesh = get_mesh(intra_elec_msh);
1194  sf_vec::ltype alg_type = sf_vec::algebraic;
1195 
1196  int dpn = 1;
1197  SF::init_vector(&kappa_i, intra_mesh, dpn, alg_type);
1198  SF::init_vector(&tmp_i1, intra_mesh, dpn, alg_type);
1199  SF::init_vector(&tmp_i2, intra_mesh, dpn, alg_type);
1200  SF::init_vector(&old_vm, intra_mesh, dpn, alg_type);
1201 
1202  if(!param_globals::operator_splitting)
1203  SF::init_vector(&Irhs, intra_mesh, dpn, alg_type);
1204 
1205  // alloc matrices
1206  int max_row_entries = max_nodal_edgecount(intra_mesh);
1207 
1211 
1212  rhs_parab->init(intra_mesh, dpn, dpn, max_row_entries);
1213  mass_i ->init(intra_mesh, dpn, dpn, param_globals::mass_lumping ? 1 : max_row_entries);
1214 
1215  parab_tech = static_cast<parabolic_solver::parabolic_t>(param_globals::parab_solve);
1216  dur = timing(t1, t0);
1217 }
1218 
1220 {
1222  double start, end, period;
1223  get_time(start);
1224  double t0, t1, dur;
1225  mass_integrator mass_integ;
1226  int dpn = 1;
1227 
1228  int log_flag = param_globals::output_level > 1 ? ECHO : 0;
1229  MaterialType & mt = mtype[Electrics::intra_grid];
1230 
1231  double Dt = user_globals::tm_manager->time_step;
1232  get_kappa(*kappa_i, param_globals::imp_region, miif, UM2_to_CM2 / Dt);
1233 
1234  cond_t condType = intra_cond;
1235  sf_mesh & mesh = get_mesh(intra_elec_msh);
1236 
1237  if( (param_globals::bidomain == MONODOMAIN && param_globals::bidm_eqv_mono) ||
1238  (param_globals::bidomain == PSEUDO_BIDM) )
1239  condType = para_cond;
1240 
1241  // set material and conductivity type
1242  set_cond_type(mt, condType);
1243 
1244  // fill the system
1245  {
1246  get_time(t0);
1247 
1248  elec_stiffness_integrator stfn_integ(mt);
1249  SF::assemble_matrix(*rhs_parab, mesh, stfn_integ);
1250 
1251  dur = timing(t1,t0);
1252  log_msg(logger,0,log_flag, "Computed parabolic stiffness matrix in %.3f seconds.", dur);
1253  get_time(t0);
1254 
1255  if(param_globals::mass_lumping)
1256  SF::assemble_lumped_matrix(*mass_i, mesh, mass_integ);
1257  else
1258  SF::assemble_matrix(*mass_i, mesh, mass_integ);
1259 
1260  sf_vec* empty; SF::init_vector(&empty);
1261  mass_i->mult_LR(*kappa_i, *empty);
1262 
1263  dur = timing(t1,t0);
1264  log_msg(logger,0,log_flag, "Computed parabolic mass matrix in %.3f seconds.", dur);
1265  delete empty;
1266  }
1267 
1268  // initialize parab lhs
1269  if(parab_tech != EXPLICIT) {
1271  // if we have mass lumping, then the nonzero pattern between Mi and Ki is different
1272  bool same_nonzero = param_globals::mass_lumping == false;
1273 
1274  if (parab_tech==CN) {
1275  lhs_parab->scale(-param_globals::theta);
1276  lhs_parab->add_scaled_matrix(*mass_i, 1.0, same_nonzero);
1277  }
1278  else if (parab_tech==O2dT) {
1279  lhs_parab->scale(-0.5);
1280  mass_i->scale(0.5);
1281  lhs_parab->add_scaled_matrix(*mass_i, 1.0, same_nonzero);
1282  lhs_parab->add_scaled_matrix(*mass_i, 1.0, same_nonzero);
1283  lhs_parab->add_scaled_matrix(*mass_i, 1.0, same_nonzero);
1284  }
1285  }
1286  else {
1289 
1290  SF_real* p = inv_mass_diag->ptr();
1291 
1292  for(int i=0; i<inv_mass_diag->lsize(); i++)
1293  p[i] = 1.0 / p[i];
1294 
1296  }
1297 
1298  if(parab_tech == CN || parab_tech == O2dT) {
1299  set_dir(INPUT);
1300  get_time(t0);
1301 
1302  setup_linear_solver(logger);
1303 
1304  dur = timing(t1,t0);
1305  log_msg(logger,0,log_flag, "Initializing parabolic solver in %.5f seconds.", dur);
1306  set_dir(OUTPUT);
1307  }
1308  period = timing(end, start);
1309 }
1310 
1311 void parabolic_solver::setup_linear_solver(FILE_SPEC logger)
1312 {
1314  tol = param_globals::cg_tol_parab;
1315  max_it = param_globals::cg_maxit_parab;
1316 
1317  std::string default_opts;
1318  std::string solver_file;
1319  solver_file = param_globals::parab_options_file;
1320  if (param_globals::flavor == std::string("ginkgo")) {
1321  default_opts = std::string(
1322  R"(
1323 {
1324  "type": "solver::Cg",
1325  "criteria": [
1326  {
1327  "type": "Iteration",
1328  "max_iters": 100
1329  },
1330  {
1331  "type": "ResidualNorm",
1332  "reduction_factor": 1e-4
1333  }
1334  ],
1335  "preconditioner": {
1336  "type": "preconditioner::Schwarz",
1337  "local_solver": {
1338  "type": "preconditioner::Ilu"
1339  }
1340  }
1341 }
1342  )");
1343  } else if (param_globals::flavor == std::string("petsc")) {
1344  default_opts = std::string("-pc_type bjacobi -sub_pc_type ilu -ksp_type cg");
1345  }
1346 
1347  lin_solver->setup_solver(*lhs_parab, tol, max_it, param_globals::cg_norm_parab,
1348  "parabolic PDE", false, logger, solver_file.c_str(),
1349  default_opts.c_str());
1350 }
1351 
1352 void parabolic_solver::solve(sf_vec & phie_i)
1353 {
1355  switch (parab_tech) {
1356  case CN: solve_CN(phie_i); break;
1357  case O2dT: solve_O2dT(phie_i); break;
1358  default: solve_EF(phie_i); break;
1359  }
1360 }
1361 
1362 void parabolic_solver::solve_CN(sf_vec & phie_i)
1363 {
1365  double t0,t1;
1366  // assembly of rhs for CN
1367  if (param_globals::bidomain == BIDOMAIN) {
1368  tmp_i1->deep_copy(phie_i);
1369  tmp_i1->add_scaled(*Vmv, 1.0 - param_globals::theta);
1370  rhs_parab->mult(*tmp_i1, *tmp_i2);
1371  }
1372  else {
1373  rhs_parab->mult(*Vmv, *tmp_i2);
1374  *tmp_i2 *= 1.0 - param_globals::theta;
1375  }
1376 
1377  mass_i->mult(*Vmv, *tmp_i1);
1378  *tmp_i1 += *tmp_i2;
1379 
1380  // add current contributions to rhs
1381  if(!param_globals::operator_splitting)
1382  tmp_i1->add_scaled(*Irhs, -1.0);
1383 
1384  get_time(t0);
1385 
1386  (*lin_solver)(*Vmv, *tmp_i1);
1387 
1388  if(lin_solver->reason < 0) {
1389  log_msg(0, 5, 0,"%s solver diverged. Reason: %s.", lin_solver->name.c_str(),
1390  petsc_get_converged_reason_str(lin_solver->reason));
1391  EXIT(1);
1392  }
1393 
1394  // treat solver statistics
1395  auto dur = timing(t1, t0);
1396  lin_solver->time += dur;
1397  stats.slvtime += dur;
1399 }
1400 
1401 void parabolic_solver::solve_O2dT(sf_vec & phie_i)
1402 {
1404  double t0,t1;
1405  // assembly of rhs for FE
1406  if (param_globals::bidomain == BIDOMAIN) {
1407  tmp_i2->deep_copy(phie_i);
1408  tmp_i2->add_scaled(*Vmv, 0.5);
1409  rhs_parab->mult(*tmp_i2, *tmp_i1); // tmp_i1 = K_i(Vm^t * 0.5 + phi_e)
1410  }
1411  else {
1412  rhs_parab->mult(*Vmv, *tmp_i1);
1413  *tmp_i1 *= 0.5; // tmp_i1 = 0.5 * K_i Vm^t
1414  }
1415 
1416  mass_i->mult(*Vmv, *tmp_i2); // tmp_i2 = M/2 Vm^t
1417  tmp_i1->add_scaled(*tmp_i2, 4.0); // tmp_i1 = (2M+K_i/2)Vm^t
1418  mass_i->mult(*old_vm, *tmp_i2); // tmp_i2 = M/2 Vm^{t-1}
1419 
1420  tmp_i1->add_scaled(*tmp_i2, -1.0); // tmp_i1 = (2M+K_i/2)Vm^t-M/2 Vm^{t-1}
1421  *old_vm = *Vmv;
1422 
1423  get_time(t0);
1424 
1425  // solve
1426  (*lin_solver)(*Vmv, *tmp_i1);
1427 
1428  // treat solver statistics
1429  stats.slvtime += timing(t1, t0);
1431 }
1432 
1433 void parabolic_solver::solve_EF(sf_vec & phie_i)
1434 {
1436  double t0,t1,t2;
1437  get_time(t0);
1438 
1439  // assembly of rhs for FE
1440  if (param_globals::bidomain == BIDOMAIN) {
1441  tmp_i2->deep_copy(phie_i);
1442  *tmp_i2 += *Vmv;
1443  rhs_parab->mult(*tmp_i2, *tmp_i1);
1444  }
1445  else {
1446  rhs_parab->mult(*Vmv, *tmp_i1);
1447  }
1448 
1449  *tmp_i1 *= *inv_mass_diag;
1450  Vmv->add_scaled(*tmp_i1, 1.0);
1451 
1452  if(param_globals::operator_splitting == false)
1453  Vmv->add_scaled(*Irhs, -1.0);
1454 
1455  // record rhs timing
1456  stats.slvtime += timing(t1, t0);
1457 }
1458 
1460 {
1461  char* prvSimDir = strlen(param_globals::start_statef) ?
1462  get_file_dir(param_globals::start_statef) : NULL;
1463 
1464  const char* extn = ".dat";
1465 
1466  // if compute_APD we need an extra 2 acts
1467  int addLATs = param_globals::compute_APD ? 2 : 0;
1468 
1469  bool have_sentinel = param_globals::t_sentinel > 0.0;
1470  bool need_to_add_sentinel = have_sentinel && (param_globals::sentinel_ID < 0);
1471 
1472  addLATs += need_to_add_sentinel ? 1 : 0;
1473  acts.resize(param_globals::num_LATs + addLATs);
1474 
1475  int j=0;
1476  for (int i = 0; i < param_globals::num_LATs; i++ )
1477  {
1478  // using Ph only with bidomain runs
1479  if (param_globals::lats[i].method <= 0 || (param_globals::lats[i].measurand == PHIE && !param_globals::bidomain)) {
1480  log_msg(NULL, 3, 0, "Phie-based LAT measurement requires bidomain >=1 Ignoring lats[%d].", i);
1481  continue;
1482  }
1483 
1484  acts[j].method = (ActMethod)param_globals::lats[i].method;
1485  acts[j].threshold = param_globals::lats[i].threshold;
1486  acts[j].mode = param_globals::lats[i].mode;
1487  acts[j].all = param_globals::lats[i].all;
1488  acts[j].measurand = (PotType)param_globals::lats[i].measurand;
1489  acts[j].ID = param_globals::lats[i].ID;
1490  acts[j].fout = NULL;
1491 
1492  if(param_globals::lats[i].all) {
1493  acts[j].fname = (char*) malloc((strlen(param_globals::lats[i].ID)+strlen(extn)+1)*sizeof(char));
1494  snprintf(acts[j].fname, strlen(param_globals::lats[i].ID)+strlen(extn)+1, "%s%s", param_globals::lats[i].ID, extn);
1495  }
1496  else {
1497  char prfx[] = "init_acts_";
1498  int max_len = strlen(prfx) + strlen(param_globals::lats[i].ID) + strlen(extn) + 1;
1499 
1500  acts[j].fname = (char*) malloc(max_len*sizeof(char));
1501  snprintf(acts[j].fname, max_len, "%s%s%s", prfx, param_globals::lats[i].ID, extn);
1502  }
1503 
1504  // restarting
1505  if(prvSimDir != NULL) {
1506  int len_fname = strlen(prvSimDir)+strlen(acts[j].fname)+2;
1507  acts[j].prv_fname = (char*) malloc(len_fname*sizeof(char));
1508  snprintf(acts[j].prv_fname, len_fname, "%s/%s", prvSimDir, acts[j].fname);
1509  }
1510 
1511  j++;
1512  }
1513 
1514  if(param_globals::compute_APD) {
1515  acts[j].method = ACT_THRESH; // threshold crossing
1516  acts[j].threshold = param_globals::actthresh;
1517  acts[j].mode = 0; // upstroke
1518  acts[j].all = true;
1519  acts[j].measurand = VM; // Vm
1520  //acts[j].ID = dupstr("Vm_Activation");
1521  acts[j].fout = NULL;
1522  acts[j].fname = dupstr("vm_activation.dat");
1523 
1524  j++;
1525  acts[j].method = ACT_THRESH; // threshold crossing
1526  acts[j].threshold = param_globals::recovery_thresh;
1527  acts[j].mode = 1; // repol
1528  acts[j].all = true;
1529  acts[j].measurand = VM; // Vm
1530  //(*acts)[j+1].ID = param_globals::lats[i].ID;
1531  acts[j].fout = NULL;
1532  acts[j].fname = dupstr("vm_repolarisation.dat");
1533 
1534  j++;
1535  }
1536 
1537  // set up sentinel for activity checking
1538  sntl.activated = have_sentinel;
1539  sntl.t_start = param_globals::t_sentinel_start;
1540  sntl.t_window = param_globals::t_sentinel;
1541  sntl.t_quiesc =-1.;
1542  sntl.ID = param_globals::sentinel_ID;
1543 
1544  if(need_to_add_sentinel) {
1545  // add a default LAT detector as sentinel
1546  acts[j].method = ACT_THRESH; // threshold crossing
1547  acts[j].threshold = param_globals::actthresh;
1548  acts[j].mode = 0; // upstroke
1549  acts[j].all = true;
1550  acts[j].measurand = VM; // Vm
1551  //(*acts)[j].ID = dupstr("Vm_Activation");
1552  acts[j].fout = NULL;
1553  acts[j].fname = dupstr("vm_sentinel.dat");
1554  // set sentinel index
1555  sntl.ID = j;
1556  j++;
1557  }
1558 
1559  if(prvSimDir) free(prvSimDir);
1560 }
1561 
1562 void print_act_log(FILE_SPEC logger, const SF::vector<Activation> & acts, int idx)
1563 {
1564  const Activation & act = acts[idx];
1565 
1566  log_msg(logger, 0, 0, "\n");
1567  log_msg(logger, 0, 0, "LAT detector [%2d]", idx);
1568  log_msg(logger, 0, 0, "-----------------\n");
1569 
1570  log_msg(logger, 0, 0, "Measurand: %s", act.measurand ? "Phie" : "Vm");
1571  log_msg(logger, 0, 0, "All: %s", act.all ? "All" : "Only first");
1572  log_msg(logger, 0, 0, "Method: %s", act.method==ACT_DT ? "Derivative" : "Threshold crossing");
1573 
1574  char buf[64], gt[2], sgn[2];
1575  snprintf(sgn, sizeof sgn, "%s", act.mode?"-":"+");
1576  snprintf(gt, sizeof gt, "%s", act.mode?"<":">");
1577 
1578  if(act.method==ACT_DT)
1579  snprintf(buf, sizeof buf, "Maximum %sdf/dt %s %.2f mV", sgn, gt, act.threshold);
1580  else
1581  snprintf(buf, sizeof buf, "Intersection %sdf/dt with %.2f", sgn, act.threshold);
1582 
1583  log_msg(logger, 0, 0, "Mode: %s", buf);
1584  log_msg(logger, 0, 0, "Threshold: %.2f mV\n", act.threshold);
1585 }
1586 
1587 void LAT_detector::init(sf_vec & vm, sf_vec & phie, int offset, enum physic_t phys_t)
1588 {
1589  if(!get_physics(phys_t)) {
1590  log_msg(0,0,5, "There seems to be no EP is defined. LAT detector requires active EP! Aborting LAT setup!");
1591  return;
1592  }
1593 
1594  // we use the electrics logger for output
1595  FILE_SPEC logger = get_physics(phys_t)->logger;
1596 
1597  // TODO(init): except for the shallow copies, shouldn't these be deleted?
1598  // When to delete them?
1599  for(size_t i = 0; i < acts.size(); ++i) {
1600  acts[i].init = 1;
1601  SF::init_vector(&(acts[i].phi));
1602  acts[i].phi->shallow_copy(!acts[i].measurand ? vm : phie);
1603  acts[i].offset = offset;
1604 
1605  SF::init_vector(&(acts[i].phip), acts[i].phi);
1606  *acts[i].phip = *acts[i].phi;
1607 
1608  // derivative based detector
1609  if (acts[i].method == ACT_DT) {
1610  SF::init_vector(&(acts[i].dvp0), acts[i].phi);
1611  SF::init_vector(&(acts[i].dvp1), acts[i].phi);
1612  if(acts[i].mode)
1613  log_msg(NULL,2,0, "Detection of -df/dt|max not implemented, +df/dt|max will be detected.");
1614  }
1615 
1616  // allocate additional local buffers
1617  acts[i].ibuf = (int *)malloc(acts[i].phi->lsize()*sizeof(int));
1618  acts[i].actbuf = (double *)malloc(acts[i].phi->lsize()*sizeof(double));
1619 
1620  if (!acts[i].all) {
1621  SF::init_vector(&acts[i].tm, acts[i].phi->gsize(), acts[i].phi->lsize());
1622  acts[i].tm->set(-1.);
1623 
1624  // initialize with previous initial activations
1625  if(acts[i].prv_fname != NULL) {
1626  set_dir(INPUT);
1627  size_t nread = acts[i].tm->read_ascii(acts[i].prv_fname);
1628  set_dir(OUTPUT);
1629 
1630  if(nread == 0) {
1631  log_msg(NULL,2,ECHO,"Warning: Initialization of LAT[%2d] failed.", i);
1632  } else {
1633  // the file is in canonical order (written with forward=true); permute back to PETSc layout
1635  assert(sc != NULL);
1636  (*sc)(*acts[i].tm, /*forward=*/false);
1637  }
1638  }
1639  }
1640  else {
1641  if ( !get_rank() ) {
1642  // here we should copy over previous file and open in append mode
1643  if(acts[i].prv_fname!=NULL) {
1644  set_dir(INPUT);
1645  FILE_SPEC in = f_open( acts[i].prv_fname, "r" );
1646  if(in) {
1647  log_msg(NULL,2,0, "Copying over of previous activation file not implemented.\n"); f_close(in);
1648  }
1649  else
1650  log_msg(NULL,3,0,"Warning: Initialization in %s - \n"
1651  "Failed to read activation file %s.\n", __func__, acts[i].prv_fname);
1652 
1653  set_dir(OUTPUT);
1654  }
1655  acts[i].fout = f_open( acts[i].fname, acts[i].prv_fname==NULL?"w":"a" );
1656  }
1657  }
1658  print_act_log(logger, acts, i);
1659  }
1660 
1661  sf_mesh & intra_mesh = get_mesh(intra_elec_msh);
1663 }
1664 
1665 
1666 int output_all_activations(FILE_SPEC fp, int *ibuf, double *act_tbuf, int nlacts)
1667 {
1668  int rank = get_rank(), gacts = 0, numProc = get_size();
1669 
1670  if (rank == 0) {
1671  // rank 0 writes directly to the table
1672  for (int i=0; i<nlacts; i++)
1673  fprintf(fp->fd, "%d\t%.6f\n", ibuf[i], act_tbuf[i]);
1674 
1675  gacts += nlacts;
1676 
1677  SF::vector<int> buf_inds;
1678  SF::vector<double> buf_acts;
1679 
1680  for (int j=1; j<numProc; j++) {
1681  int acts = 0;
1682  MPI_Status status;
1683  MPI_Recv(&acts, 1, MPI_INT, j, 110, PETSC_COMM_WORLD, &status);
1684 
1685  if (acts) {
1686  buf_inds.resize(acts);
1687  buf_acts.resize(acts);
1688 
1689  MPI_Recv(buf_inds.data(), acts, MPI_INT, j, 110, PETSC_COMM_WORLD, &status);
1690  MPI_Recv(buf_acts.data(), acts, MPI_DOUBLE, j, 110, PETSC_COMM_WORLD, &status);
1691 
1692  for(int ii=0; ii<acts; ii++)
1693  fprintf(fp->fd, "%d\t%.6f\n", buf_inds[ii], buf_acts[ii]);
1694 
1695  gacts += acts;
1696  }
1697  }
1698  fflush(fp->fd);
1699  }
1700  else {
1701  MPI_Send(&nlacts, 1, MPI_INT, 0, 110, PETSC_COMM_WORLD);
1702  if (nlacts) {
1703  MPI_Send(ibuf, nlacts, MPI_INT, 0, 110, PETSC_COMM_WORLD);
1704  MPI_Send(act_tbuf, nlacts, MPI_DOUBLE, 0, 110, PETSC_COMM_WORLD);
1705  }
1706  }
1707 
1708  MPI_Bcast(&gacts, 1, MPI_INT, 0, PETSC_COMM_WORLD);
1709  return gacts;
1710 }
1712 int LAT_detector::check_acts(double tm)
1713 {
1714  int nacts = 0;
1715  double *a;
1716 
1717  for(Activation* aptr = acts.data(); aptr != acts.end(); aptr++)
1718  {
1719  int lacts = 0;
1720  switch (aptr->method) {
1721  case ACT_THRESH:
1722  lacts = check_cross_threshold(*aptr->phi, *aptr->phip, tm,
1723  aptr->ibuf, aptr->actbuf, aptr->threshold, aptr->mode);
1724  break;
1725 
1726  case ACT_DT:
1727  lacts = check_mx_derivative (*aptr->phi, *aptr->phip, tm,
1728  aptr->ibuf, aptr->actbuf, *aptr->dvp0, *aptr->dvp1,
1729  aptr->threshold, aptr->mode);
1730  break;
1731 
1732  default:
1733  break;
1734  }
1735 
1736  if (!aptr->all)
1737  a = aptr->tm->ptr();
1738 
1740 
1741  for(int j=0; j<lacts; j++) {
1742  if(aptr->all) {
1743  int nodal_idx = this->petsc_to_nodal.forward_map(aptr->ibuf[j]);
1744  aptr->ibuf[j] = canon_nbr[nodal_idx] + aptr->offset;
1745  }
1746  else {
1747  if(a[aptr->ibuf[j]] == -1)
1748  a[aptr->ibuf[j]] = aptr->actbuf[j];
1749  }
1750  }
1751 
1752  if(aptr->all)
1753  output_all_activations(aptr->fout, aptr->ibuf, aptr->actbuf, lacts);
1754  else
1755  aptr->tm->release_ptr(a);
1756 
1757  MPI_Allreduce(MPI_IN_PLACE, &lacts, 1, MPI_INT, MPI_SUM, PETSC_COMM_WORLD);
1758  nacts += lacts;
1759 
1760  aptr->nacts = nacts;
1761  }
1762 
1763  return nacts > 0;
1764 }
1765 
1766 
1767 int LAT_detector::check_quiescence(double tm, double dt)
1768 {
1769  static int savequitFlag = 0;
1770  int numNodesActivated = -1;
1771 
1772  if(sntl.activated) {
1773  // initialization
1774  if(sntl.t_quiesc < 0. && sntl.t_window >= 0.0 ) {
1775  log_msg(0,0,ECHO | NONL, "================================================================================================\n");
1776  log_msg(0,0,ECHO | NONL, "%s() WARNING: simulation is configured to savequit() after %.2f ms of quiescence\n", __func__, sntl.t_window);
1777  log_msg(0,0,ECHO | NONL, "================================================================================================\n");
1778  sntl.t_quiesc = 0.0;
1779  }
1780 
1781  if(tm >= sntl.t_start && !savequitFlag)
1782  {
1783  numNodesActivated = acts[sntl.ID].nacts;
1784 
1785  if(numNodesActivated) sntl.t_quiesc = 0.0;
1786  else sntl.t_quiesc += dt;
1787 
1788  if(sntl.t_window >= 0.0 && sntl.t_quiesc > sntl.t_window && !savequitFlag) {
1789  savequitFlag++;
1790  savequit();
1791  }
1792  }
1793  }
1794 
1795  return numNodesActivated;
1796 }
1797 
1798 
1799 
1800 
1801 int LAT_detector::check_cross_threshold(sf_vec & vm, sf_vec & vmp, double tm,
1802  int *ibuf, double *actbuf, float threshold, int mode)
1803 {
1804  SF_real *c = vm.ptr();
1805  SF_real *p = vmp.ptr();
1806  int lsize = vm.lsize();
1807  int nacts = 0, gnacts = 0;
1808 
1809  for (int i=0; i<lsize; i++) {
1810  int sgn = 1;
1811  bool triggered = false;
1812  if(mode==0) {// detect +slope crossing
1813  triggered = p[i] <= threshold && c[i] > threshold; }
1814  else { // detect -slope crossing
1815  triggered = p[i] >= threshold && c[i] < threshold;
1816  sgn = -1;
1817  }
1818 
1819  if (triggered) {
1820  double tact = tm - param_globals::dt + (threshold-p[i])/(c[i]-p[i])*sgn*param_globals::dt;
1821  ibuf [nacts] = i;
1822  actbuf[nacts] = tact;
1823  nacts++;
1824  }
1825  p[i] = c[i];
1826  }
1827 
1828  vm.release_ptr(c);
1829  vmp.release_ptr(p);
1830  return nacts;
1831 }
1832 
1833 int LAT_detector::check_mx_derivative(sf_vec & vm, sf_vec & vmp, double tm,
1834  int *ibuf, double *actbuf, sf_vec & dvp0, sf_vec & dvp1,
1835  float threshold, int mode)
1836 {
1837  int nacts = 0, gnacts = 0;
1838  double tact, dt2 = 2 * param_globals::dt;
1839  int lsize = vm.lsize();
1840  SF_real ddv0, ddv1, dv, dvdt;
1841  SF_real *c, *p, *pd0, *pd1;
1842 
1843  c = vm.ptr();
1844  p = vmp.ptr();
1845  pd0 = dvp0.ptr();
1846  pd1 = dvp1.ptr();
1847 
1848  for (int i=0; i<lsize; i++ ) {
1849  dv = (c[i]-p[i]);
1850  dvdt = dv/param_globals::dt;
1851  ddv0 = pd1[i]-pd0[i];
1852  ddv1 = dv -pd1[i];
1853 
1854  if (dvdt>=threshold && ddv0>0 && ddv1<0) {
1855  tact = tm-dt2+(ddv0/(ddv0-ddv1))*param_globals::dt;
1856  ibuf [nacts] = i;
1857  actbuf[nacts] = tact;
1858  nacts++;
1859  }
1860  p [i] = c[i];
1861  pd0[i] = pd1[i];
1862  pd1[i] = dv;
1863  }
1864 
1865  vm .release_ptr(c);
1866  vmp .release_ptr(p);
1867  dvp0.release_ptr(pd0);
1868  dvp1.release_ptr(pd1);
1869 
1870  return nacts;
1871 }
1872 
1877 {
1879  assert(sc != NULL);
1880 
1881  bool forward = true;
1882 
1883  for (size_t i = 0; i < acts.size(); i++) {
1884  if (is_init(acts[i].tm)) {
1885  (*sc)(*acts[i].tm, forward);
1886  acts[i].tm->write_ascii(acts[i].fname, false);
1887  }
1888  }
1889 }
1890 
1891 void Electrics::prepace() {
1892  log_msg(NULL, 0, 0, "Using activation times from file %s to distribute prepacing states\n",
1893  param_globals::prepacing_lats);
1894  log_msg(NULL, 0, 0, "Assuming stimulus strength %f uA/uF with duration %f ms for prepacing\n",
1895  param_globals::prepacing_stimstr, param_globals::prepacing_stimdur);
1896 
1897  limpet::MULTI_IF* miif = this->ion.miif;
1898 
1899  const sf_mesh & mesh = get_mesh(intra_elec_msh);
1900  sf_vec* read_lats; SF::init_vector(&read_lats, mesh, 1, sf_vec::algebraic);
1901 
1902  // read in the global distributed vector of all activation times
1903  set_dir(INPUT);
1904  size_t numread = read_lats->read_ascii(param_globals::prepacing_lats);
1905  if (numread == 0) {
1906  log_msg(NULL, 5, 0, "Failed reading required LATs! Skipping prepacing!");
1907  return;
1908  }
1909  set_dir(OUTPUT);
1910 
1912  assert(sc != NULL);
1913 
1914  // permute in-place to petsc permutation
1915  bool forward = false;
1916  (*sc)(*read_lats, forward);
1917 
1918  // take care of negative LAT values
1919  {
1920  PetscReal* lp = read_lats->ptr();
1921  for(int i=0; i<read_lats->lsize(); i++)
1922  if(lp[i] < 0.0) lp[i] = param_globals::tend + 10.0;
1923 
1924  read_lats->release_ptr(lp);
1925  }
1926 
1927  // make LATs relative and figure out the first LAT
1928  // so we know when to save state of each point
1929  SF_real LATmin = read_lats->min();
1930 
1931  if(LATmin < 0.0) {
1932  log_msg(0,3,0, "LAT data is not complete. Skipping prepacing.");
1933  return;
1934  }
1936  SF_real offset = floor(LATmin / param_globals::prepacing_bcl) * param_globals::prepacing_bcl;
1937  SF_real last_tm = param_globals::prepacing_bcl * param_globals::prepacing_beats;
1938 
1939  // compute read_lats[i] = last_tm - (read_lats[i] - offset)
1940  *read_lats += -offset;
1941  *read_lats *= -1.;
1942  *read_lats += last_tm;
1943 
1944  miif->getRealData();
1945  SF_real *save_tm = read_lats->ptr();
1946  SF_real *vm = miif->gdata[limpet::Vm]->ptr();
1947 
1948  for (int ii = 0; ii < miif->N_IIF; ii++) {
1949  if (!miif->N_Nodes[ii]) continue;
1950 
1951  // create sorted array of save times.
1952  SF::vector<SF::mixed_tuple<double,int>> sorted_save(miif->N_Nodes[ii]); // v1 = time, v2 = index
1953  for (int kk = 0; kk < miif->N_Nodes[ii]; kk++) {
1954  sorted_save[kk].v1 = save_tm[miif->NodeLists[ii][kk]];
1955  sorted_save[kk].v2 = kk;
1956  }
1957  std::sort(sorted_save.begin(), sorted_save.end());
1958 
1959  size_t lastidx = sorted_save.size() - 1;
1960  int paced = sorted_save[lastidx].v2; // IMP index of latest node
1961  int csav = 0;
1962 
1963  for (double t = 0; t < sorted_save[lastidx].v1; t += param_globals::dt) {
1964  if (fmod(t, param_globals::prepacing_bcl) < param_globals::prepacing_stimdur &&
1965  t < param_globals::prepacing_bcl * param_globals::prepacing_beats - 1)
1966  miif->ldata[ii][limpet::Vm][paced] += param_globals::prepacing_stimstr * param_globals::dt;
1967 
1968  compute_IIF(*miif->IIF[ii], miif->ldata[ii], paced);
1969 
1970  // Vm update always happens now outside of the imp
1971  miif->ldata[ii][limpet::Vm][paced] -= miif->ldata[ii][limpet::Iion][paced] * param_globals::dt;
1972  vm[miif->NodeLists[ii][paced]] = miif->ldata[ii][limpet::Vm][paced];
1973 
1974  while (csav < miif->N_Nodes[ii] - 1 && t >= sorted_save[csav].v1)
1975  dup_IMP_node_state(*miif->IIF[ii], paced, sorted_save[csav++].v2, miif->ldata[ii]);
1976  }
1977 
1978  // get nodes which may be tied for last
1979  while (csav < miif->N_Nodes[ii] - 1)
1980  dup_IMP_node_state(*miif->IIF[ii], paced, sorted_save[csav++].v2, miif->ldata[ii]);
1981  // ipdate global Vm vector
1982  for (int k = 0; k < miif->N_Nodes[ii]; k++) vm[miif->NodeLists[ii][k]] = miif->ldata[ii][limpet::Vm][k];
1983  }
1984 
1985  read_lats->release_ptr(save_tm);
1986  miif->gdata[limpet::Vm]->release_ptr(vm);
1987  miif->releaseRealData();
1988 }
1989 
1990 
1991 void recover_phie_std(sf_vec & vm, phie_recovery_data & rcv)
1992 {
1994  if (!rcv.pts.size())
1995  return;
1996 
1997  int rank = get_rank();
1998 
1999  if(!get_physics(elec_phys)) {
2000  log_msg(0,0,5, "There seems to be no EP is defined. Phie recovery requires active EP! Aborting!");
2001  return;
2002  }
2003 
2004  Electrics* elec = static_cast<Electrics*>(get_physics(elec_phys));
2005  sf_mat & Ki = *elec->parab_solver.rhs_parab;
2006 
2007  const sf_mesh & imesh = get_mesh(intra_elec_msh);
2008  const SF::vector<mesh_int_t> & alg_nod = imesh.pl.algebraic_nodes();
2009 
2010  SF_int start, end;
2011  vm.get_ownership_range(start, end);
2012 
2013  if(!rcv.Im) {
2014  SF::init_vector(&rcv.Im, &vm);
2015  SF::init_vector(&rcv.dphi, &vm);
2016  }
2017 
2018  SF_int r_start, r_end;
2019  rcv.phie_rec->get_ownership_range(r_start, r_end);
2020 
2021  SF_real *ph_r = rcv.phie_rec->ptr();
2022 
2023  // use minimum distance to ensure r>0
2024  // consistent with the line source approximation, the "cable radius"
2025  // is used as a lower limit for the source-field point distance
2026  float minDist = 2. / param_globals::imp_region[0].cellSurfVolRatio; // radius in um
2027 
2028  Ki.mult(vm, *rcv.Im);
2029  int numpts = rcv.pts.size() / 3;
2030  Point fpt, cpt;
2031 
2032  for (int j=0; j<numpts; j++) {
2033  fpt = rcv.pts.data() + j*3;
2034 
2035  *rcv.dphi = *rcv.Im;
2036  SF_real* dp = rcv.dphi->ptr();
2037 
2038  for (size_t i = 0; i<alg_nod.size(); i++)
2039  {
2040  mesh_int_t loc_nodal_idx = alg_nod[i];
2041  mesh_int_t loc_petsc_idx = local_nodal_to_local_petsc(imesh, rank, loc_nodal_idx);
2042  cpt = imesh.xyz.data()+loc_nodal_idx*3;
2043 
2044  double r = dist(fpt, cpt) + minDist;
2045  dp[loc_petsc_idx] /= r;
2046  }
2047 
2048  rcv.dphi->release_ptr(dp);
2049 
2050  SF_real phi = rcv.dphi->sum() / 4. / M_PI / rcv.gBath;
2051  if ( (j>=r_start) && (j<r_end) )
2052  ph_r[j-r_start] = phi;
2053  }
2054 
2055  rcv.phie_rec->release_ptr(ph_r);
2056 }
2057 
2059 {
2060  int err = 0, rank = get_rank();
2061 
2063  log_msg(0,0,5, "There seems to be no EP is defined. Phie recovery requires active EP! Aborting!");
2064  return 1;
2065  }
2066 
2067  sf_mesh & imesh = get_mesh(intra_elec_msh);
2068  Electrics* elec = static_cast<Electrics*>(get_physics(elec_phys));
2069  phie_recovery_data & phie_rcv = elec->phie_rcv;
2070 
2071  // we close the files of the default electrics if there are any open
2072  elec->output_manager.close_files_and_cleanup();
2073 
2074  // register output
2075  set_dir(POSTPROC);
2076  igb_output_manager phie_rec_out;
2077  phie_rec_out.register_output(phie_rcv.phie_rec, phie_recv_msh, 1,
2078  param_globals::phie_recovery_file, "mV");
2079 
2080  // Buffer for Vm data
2081  sf_vec* vm = get_data(vm_vec); assert(vm);
2082 
2083  // set up igb header and point fd to start of Vm file
2084  set_dir(OUTPUT);
2085  IGBheader vm_igb;
2086  if(rank == 0) {
2087  FILE_SPEC file = f_open(param_globals::vofile, "r");
2088  if(file != NULL) {
2089  vm_igb.fileptr(file->fd);
2090  vm_igb.read();
2091 
2092  if(vm_igb.x() != vm->gsize()) {
2093  log_msg(0,4,0, "%s error: Vm dimension does not fit to %s file. Aborting recovery! \n",
2094  __func__, param_globals::vofile);
2095  err++;
2096  }
2097 
2098  delete file;
2099  }
2100  else err++;
2101  }
2102 
2103  err = get_global(err, MPI_MAX);
2104 
2105  if(err == 0) {
2106  FILE* fd = static_cast<FILE*>(vm_igb.fileptr());
2107 
2108  // number of data slices
2109  const int num_io = user_globals::tm_manager->timers[iotm_spacedt]->numIOs;
2110 
2111  // scatterers
2113  assert(petsc_to_canonical != NULL);
2114 
2115  // loop over vm slices and recover phie
2116  for(int i=0; i<num_io; i++) {
2117  log_msg(0,0,0, "Step %d / %d", i+1, num_io);
2118  size_t nread = vm->read_binary<float>(fd);
2119 
2120  if(nread != size_t(vm->gsize())) {
2121  log_msg(0,3,0, "%s warning: read incomplete data slice! Aborting!", __func__);
2122  err++;
2123  break;
2124  }
2125 
2126  // permute vm_buff
2127  bool forward = false;
2128  (*petsc_to_canonical)(*vm, forward);
2129 
2130  // do phie computation
2131  recover_phie_std(*vm, phie_rcv);
2132 
2133  phie_rec_out.write_data();
2134  }
2135 
2136  phie_rec_out.close_files_and_cleanup();
2137  }
2138  return err;
2139 }
2140 
2141 void setup_phie_recovery_data(phie_recovery_data & data)
2142 {
2144  if(!get_physics(elec_phys) ) {
2145  log_msg(0,0,5, "There seems to be no EP is defined. Phie recovery requires active EP! Aborting!");
2146  return;
2147  }
2148 
2149  int rank = get_rank(), size = get_size();
2150  Electrics* elec = static_cast<Electrics*>(get_physics(elec_phys));
2151 
2152  sf_mesh & imesh = get_mesh(intra_elec_msh);
2153  const std::string basename = param_globals::phie_rec_ptf;
2154  SF::vector<mesh_int_t> ptsidx;
2155 
2156  set_dir(INPUT);
2157  SF::read_points(basename, imesh.comm, data.pts, ptsidx);
2158  make_global(data.pts, imesh.comm); // we want all ranks to have all points
2159 
2160  // set up parallel layout of recovery points
2161  SF::vector<mesh_int_t> layout;
2162  layout_from_count(mesh_int_t(ptsidx.size()), layout, imesh.comm);
2163 
2164  // set up petsc_vector for recovered potentials
2165  SF::init_vector(&data.phie_rec, layout[size], layout[rank+1]-layout[rank], 1, sf_vec::algebraic);
2166 
2167  // get conductivty
2168  SF::vector<RegionSpecs> & intra_regions = elec->mtype[Electrics::intra_grid].regions;
2169  data.gBath = static_cast<elecMaterial*>(intra_regions[0].material)->BathVal[0];
2170 }
2171 
2172 void Laplace::initialize()
2173 {
2175  int rank = get_rank();
2176 
2177  assert(param_globals::bidomain == BIDOMAIN);
2178  double t1, t2;
2179  get_time(t1);
2180 
2181  // set up Extracellular tissue
2184  mtype[Electrics::extra_grid].regionIDs, true, "gregion_e");
2185 
2186  // set up a subset of the complete electrical mappings
2187  int dpn = 1;
2189 
2191  // set up Intracellular tissue
2194  mtype[Electrics::intra_grid].regionIDs, true, "gregion_i");
2195 
2198  }
2199 
2200  // set up stimuli
2201  init_stim_info();
2202  stimuli.resize(param_globals::num_stim);
2204  for(int i=0; i<param_globals::num_stim; i++) {
2205  // construct new stimulus
2206  stimulus & s = stimuli[i];
2207 
2209  s.translate(i);
2210 
2211  s.setup(i);
2212 
2213  if(s.phys.type == Phi_ex) {
2214  s.pulse.wform = constPulse;
2215  sample_wave_form(s.pulse, i);
2216  }
2217  }
2218 
2219  set_dir(OUTPUT);
2220 
2221  ellip_solver.init();
2223 
2224  if(param_globals::dump2MatLab) {
2225  std::string bsname = param_globals::dump_basename;
2226  std::string fn;
2227 
2228  set_dir(OUTPUT);
2229  fn = bsname + "_Kie.bin";
2230  ellip_solver.phie_mat->write(fn.c_str());
2231  }
2232 
2233  // the laplace solver executes only once, thus we need a singlestep timer
2234  timer_idx = user_globals::tm_manager->add_singlestep_timer(0.0, 0.0, "laplace trigger", nullptr);
2235 
2236  SF::vector<mesh_int_t>* restr_i = NULL;
2237  SF::vector<mesh_int_t>* restr_e = NULL;
2238 
2239  setup_dataout(param_globals::dataout_e, param_globals::dataout_e_vtx, extra_elec_msh,
2240  restr_e, param_globals::num_io_nodes > 0);
2241  if(param_globals::dataout_e)
2242  output_manager.register_output(ellip_solver.phie, extra_elec_msh, 1, param_globals::phiefile, "mV", restr_e);
2243 
2245  setup_dataout(param_globals::dataout_i, param_globals::dataout_i_vtx, intra_elec_msh,
2246  restr_i, param_globals::num_io_nodes > 0);
2247  if(param_globals::dataout_i)
2248  output_manager.register_output(ellip_solver.phie_i, intra_elec_msh, 1, param_globals::phieifile, "mV", restr_i);
2249  }
2250 
2251  this->initialize_time += timing(t2, t1);
2253  this->compute_step();
2254 }
2255 
2256 void Laplace::destroy()
2257 {}
2258 
2259 void Laplace::compute_step()
2260 {
2262  // Laplace compute might be called multiple times, we want to run only once..
2263  if(!ellip_solver.lin_solver) return;
2264 
2265  double t0, t1, dur;
2266  log_msg(0,0,0, "Solving Laplace problem ..");
2267 
2268  get_time(t0);
2270  dur = timing(t1,t0);
2271 
2272  log_msg(0,0,0, "Done in %.5f seconds.", dur);
2273 
2275  this->compute_time += timing(t1, t0);
2276  set_dir(OUTPUT);
2279 
2280  // we clear the elliptic matrices and solver to save some memory when computing
2281  // the laplace solution on-the-fly
2282  delete ellip_solver.mass_e; ellip_solver.mass_e = NULL;
2283  delete ellip_solver.phie_mat; ellip_solver.phie_mat = NULL;
2285 }
2286 
2287 void Laplace::output_step()
2288 {}
2289 
2290 double Laplace::timer_val(const int timer_id)
2291 {
2292  int sidx = stimidx_from_timeridx(stimuli, timer_id);
2293  double val = 0.0;
2294 
2295  if(sidx != -1) stimuli[sidx].value(val);
2296  else val = std::nan("NaN");
2297  return val;
2298 }
2299 
2300 std::string Laplace::timer_unit(const int timer_id)
2301 {
2302  int sidx = stimidx_from_timeridx(stimuli, timer_id);
2303  std::string s_unit;
2304  if(sidx != -1) s_unit = stimuli[sidx].pulse.wave.f_unit;
2305  return s_unit;
2306 }
2307 
2309  sf_mat & mass_i,
2310  sf_mat & mass_e,
2311  limpet::MULTI_IF *miif,
2312  FILE_SPEC logger)
2313 {
2315 
2316  for(stimulus & s : stimuli) {
2317  if(is_current(s.phys.type) && s.phys.total_current) {
2318  // extracellular current injection
2319  if (s.phys.type == I_ex) {
2320  // compute affected volume in um^3
2321  SF_real vol = get_volume_from_nodes(mass_e, s.electrode.vertices);
2322 
2323  // s->strength holds the total current in uA, compute current density in uA/cm^3
2324  // Theoretically, we don't need to scale the volume to cm^3 here since we later
2325  // multiply with the mass matrix and we get um^3 * uA/um^3 = uA.
2326  // However, for I_ex there is an additional um^3 to cm^3 scaling in phys.scale,
2327  // since I_e is expected to be in uA/cm^3. Therefore, we need to compensate for that to arrive at uA later.
2328  float scale = 1.e12/vol;
2329 
2330  s.pulse.strength *= scale;
2331 
2332  log_msg(logger,0,ECHO,
2333  "%s [Stimulus %d]: current density scaled to %.4g uA/cm^3\n",
2334  s.name.c_str(), s.idx, s.pulse.strength);
2335  }
2336  else if (s.phys.type == I_tm) {
2337  // compute affected volume in um^3
2338  SF_real vol = get_volume_from_nodes(mass_i, s.electrode.vertices);
2339  const sf_mesh & imesh = get_mesh(intra_elec_msh);
2340  const SF::vector<mesh_int_t> & alg_nod = imesh.pl.algebraic_nodes();
2341 
2342  if(alg_idx_map.size() == 0) {
2343  mesh_int_t lidx = 0;
2344  for(mesh_int_t n : alg_nod) {
2345  alg_idx_map[n] = lidx;
2346  lidx++;
2347  }
2348  }
2349 
2350  SF_real surf = 0.0;
2351  for(mesh_int_t n : s.electrode.vertices) {
2352  if(alg_idx_map.count(n)) {
2353  mesh_int_t lidx = alg_idx_map[n];
2354  int r = miif->IIFmask[lidx];
2355  // surf = vol*beta [1/um], surf is in [um^2]
2356  surf = vol * miif->IIF[r]->cgeom().SVratio * param_globals::imp_region[r].volFrac;
2357  //convert to cm^2
2358  surf /= 1.e8;
2359  break;
2360  }
2361  }
2362  surf = get_global(surf, MPI_MAX, PETSC_COMM_WORLD);
2363 
2364  // scale surface density now to result in correct total current
2365  s.pulse.strength /= surf;
2366  log_msg(logger, 0, ECHO,
2367  "%s [Stimulus %d]: current density scaled to %.4g uA/cm^2\n",
2368  s.name.c_str(), s.idx, s.pulse.strength);
2369  }
2370  }
2371  }
2372 }
2373 
2374 
2375 
2376 } // namespace opencarp
#define M_PI
Definition: ION_IF.h:52
int mesh_int_t
Definition: SF_container.h:46
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
std::int32_t SF_int
Use the general std::int32_t as int type.
Definition: SF_globals.h:37
#define ECHO
Definition: basics.h:324
#define NONL
Definition: basics.h:328
#define CALI_CXX_MARK_FUNCTION
Definition: caliper_hooks.h:5
virtual void mult(const abstract_vector< T, S > &x, abstract_vector< T, S > &b) const =0
virtual void scale(S s)=0
virtual void zero()=0
virtual void get_diagonal(abstract_vector< T, S > &vec) const =0
virtual void mult_LR(const abstract_vector< T, S > &L, const abstract_vector< T, S > &R)=0
virtual void init(T iNRows, T iNCols, T ilrows, T ilcols, T loc_offset, T mxent)
virtual void duplicate(const abstract_matrix< T, S > &M)=0
virtual void add_scaled_matrix(const abstract_matrix< T, S > &A, const S s, const bool same_nnz)=0
virtual void write(const char *filename) const =0
virtual S mag() const =0
virtual S * ptr()=0
virtual void release_ptr(S *&p)=0
virtual void deep_copy(const abstract_vector< T, S > &v)=0
virtual void shallow_copy(const abstract_vector< T, S > &v)=0
virtual void add_scaled(const abstract_vector< T, S > &vec, S k)=0
virtual T lsize() const =0
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false, const bool local=false)=0
const meshdata< mesh_int_t, mesh_real_t > * mesh
the connected mesh
T forward_map(T idx) const
Map one index from a to b.
Definition: SF_container.h:264
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:429
vector< T > con
Definition: SF_container.h:412
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:464
Container for a PETSc VecScatter.
void forward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Forward scattering.
void backward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Backward scattering.
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
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
T & push_back(T val)
Definition: SF_vector.h:283
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:623
size_t size() const
Definition: hashmap.hpp:731
int numNode
local number of nodes
Definition: MULTI_ION_IF.h:210
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 dump_state(char *, float, opencarp::mesh_t gid, bool, unsigned int)
GlobalData_t *** ldata
data local to each IMP
Definition: MULTI_ION_IF.h:205
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
IIF_Mask_t * IIFmask
region for each node
Definition: MULTI_ION_IF.h:214
int timer_idx
the timer index received from the timer manager
Definition: physics_types.h:66
FILE_SPEC logger
The logger of the physic, each physic should have one.
Definition: physics_types.h:64
SF::vector< stimulus > stimuli
the electrical stimuli
Definition: electrics.h:263
LAT_detector lat
the activation time detector
Definition: electrics.h:275
grid_t
An electrics grid identifier to distinguish between intra and extra grids.
Definition: electrics.h:257
phie_recovery_data phie_rcv
struct holding helper data for phie recovery
Definition: electrics.h:284
generic_timing_stats IO_stats
Definition: electrics.h:286
void destroy()
Currently we only need to close the file logger.
Definition: electrics.cc:370
gvec_data gvec
datastruct holding global IMP state variable output
Definition: electrics.h:278
elliptic_solver ellip_solver
Solver for the elliptic bidomain equation.
Definition: electrics.h:270
MaterialType mtype[2]
the material types of intra_grid and extra_grid grids.
Definition: electrics.h:261
std::string timer_unit(const int timer_id)
figure out units of a signal linked to a given timer
Definition: electrics.cc:767
parabolic_solver parab_solver
Solver for the parabolic bidomain equation.
Definition: electrics.h:272
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
Definition: electrics.cc:751
void initialize()
Initialize the Electrics.
Definition: electrics.cc:43
igb_output_manager output_manager
class handling the igb output
Definition: electrics.h:281
int read(bool quiet=false)
Definition: IGBheader.cc:761
void fileptr(gzFile f)
Definition: IGBheader.cc:336
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
SF::index_mapping< mesh_int_t > petsc_to_nodal
Definition: electrics.h:220
int check_quiescence(double tm, double dt)
check for quiescence
Definition: electrics.cc:1711
void output_initial_activations()
output one nodal vector of initial activation time
Definition: electrics.cc:1820
void init(sf_vec &vm, sf_vec &phie, int offset, enum physic_t=elec_phys)
initializes all datastructs after electric solver setup
Definition: electrics.cc:1531
int check_acts(double tm)
check activations at sim time tm
Definition: electrics.cc:1656
SF::vector< Activation > acts
Definition: electrics.h:219
LAT_detector()
constructor, sets up basic datastructs from global_params
Definition: electrics.cc:1403
SF::vector< stimulus > stimuli
the electrical stimuli
Definition: electrics.h:390
elliptic_solver ellip_solver
Solver for the elliptic bidomain equation.
Definition: electrics.h:393
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
Definition: electrics.cc:2234
std::string timer_unit(const int timer_id)
figure out units of a signal linked to a given timer
Definition: electrics.cc:2244
MaterialType mtype[2]
the material types of intra_grid and extra_grid grids.
Definition: electrics.h:388
igb_output_manager output_manager
class handling the igb output
Definition: electrics.h:395
manager for dirichlet boundary conditions
Definition: stimulate.h:205
void enforce_dbc_rhs(sf_vec &rhs)
Definition: stimulate.cc:675
void recompute_dbcs()
recompute the dbc data.
Definition: stimulate.cc:596
bool dbc_update()
check if dbcs have updated
Definition: stimulate.cc:635
sf_mat * phie_mat
lhs matrix to solve elliptic
Definition: electrics.h:54
void rebuild_stiffness(MaterialType *mtype, SF::vector< stimulus > &stimuli, FILE_SPEC logger)
Definition: electrics.cc:953
void rebuild_matrices(MaterialType *mtype, SF::vector< stimulus > &stimuli, FILE_SPEC logger)
Definition: electrics.cc:941
lin_solver_stats stats
Definition: electrics.h:60
sf_vec * phie_i
phi_e on intracellular grid
Definition: electrics.h:49
void solve(sf_mat &Ki, sf_vec &Vmv, sf_vec &tmp_i)
Definition: electrics.cc:1058
sf_vec * phie
phi_e
Definition: electrics.h:48
sf_sol * lin_solver
petsc or ginkgo lin_solver
Definition: electrics.h:57
sf_mat * mass_e
mass matrix for RHS elliptic calc
Definition: electrics.h:53
double tol
CG stopping tolerance.
Definition: electrics.h:67
sf_vec * currtmp
temp vector for phiesrc
Definition: electrics.h:51
dbc_manager * dbc
dbcs require a dbc manager
Definition: electrics.h:63
int max_it
maximum number of iterations
Definition: electrics.h:68
sf_vec * phiesrc
I_e.
Definition: electrics.h:50
void rebuild_mass(FILE_SPEC logger)
Definition: electrics.cc:1014
void write_data()
write registered data to disk
Definition: sim_utils.cc:1537
void register_output_sync(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)
Definition: sim_utils.cc:1373
void close_files_and_cleanup()
close file descriptors
Definition: sim_utils.cc:1584
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:1504
sf_vec * Ivol
global Vm vector
Definition: electrics.h:107
double tol
CG stopping tolerance.
Definition: electrics.h:132
sf_vec * Iact
global Vm vector
Definition: electrics.h:108
sf_vec * Diff_term
Diffusion current.
Definition: electrics.h:123
sf_mat * rhs_parab
rhs matrix to solve parabolic
Definition: electrics.h:119
sf_vec * kappa_i
scaling vector for intracellular mass matrix, M
Definition: electrics.h:111
lin_solver_stats stats
Definition: electrics.h:129
void rebuild_matrices(MaterialType *mtype, limpet::MULTI_IF &miif, FILE_SPEC logger)
Definition: electrics.cc:1183
parabolic_t parab_tech
manner in which parabolic equations are solved
Definition: electrics.h:134
void solve(sf_vec &phie_i)
Definition: electrics.cc:1296
sf_vec * inv_mass_diag
inverse diagonal of mass matrix, for EXPLICIT solving
Definition: electrics.h:115
sf_mat * mass_i
lumped for parabolic problem
Definition: electrics.h:118
sf_vec * Ic
global Vm vector
Definition: electrics.h:106
sf_vec * tmp_i2
scratch vector for i-grid
Definition: electrics.h:113
int max_it
maximum number of iterations
Definition: electrics.h:133
sf_vec * tmp_i1
scratch vector for i-grid
Definition: electrics.h:112
sf_mat * lhs_parab
lhs matrix (CN) to solve parabolic
Definition: electrics.h:120
sf_vec * Vmv
global Vm vector
Definition: electrics.h:104
sf_vec * Irhs
weighted transmembrane currents
Definition: electrics.h:114
sf_vec * old_vm
older Vm needed for 2nd order dT
Definition: electrics.h:110
sf_sol * lin_solver
petsc or ginkgo lin_solver
Definition: electrics.h:126
sf_vec * IIon
ionic currents
Definition: electrics.h:103
SF::vector< mesh_int_t > vertices
Definition: stimulate.h:152
bool total_current
whether we apply total current scaling
Definition: stimulate.h:140
stim_t type
type of stimulus
Definition: stimulate.h:137
int timer_id
timer for stimulus
Definition: stimulate.h:122
waveform_t wform
wave form of stimulus
Definition: stimulate.h:95
double strength
strength of stimulus
Definition: stimulate.h:93
stim_protocol ptcl
applied stimulation protocol used
Definition: stimulate.h:168
int idx
index in global input stimulus array
Definition: stimulate.h:164
stim_electrode electrode
electrode geometry
Definition: stimulate.h:170
stim_pulse pulse
stimulus wave form
Definition: stimulate.h:167
void translate(int id)
convert legacy definitions to new format
Definition: stimulate.cc:105
void setup(int idx)
Setup from a param stimulus index.
Definition: stimulate.cc:166
stim_physics phys
physics of stimulus
Definition: stimulate.h:169
bool value(double &v) const
Get the current value if the stimulus is active.
Definition: stimulate.cc:437
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
int add_singlestep_timer(double tg, double idur, const char *iname, const char *poolname=nullptr)
Definition: timer_utils.h:143
long d_end
final index in multiples of dt
Definition: timer_utils.h:82
std::vector< base_timer * > timers
vector containing individual timers
Definition: timer_utils.h:84
double time
current time
Definition: timer_utils.h:76
Tissue level electrics, main Electrics physics class.
#define DUMP_IC
Definition: electrics.h:39
#define DUMP_IACT
Definition: electrics.h:41
#define DUMP_IVOL
Definition: electrics.h:40
void init_solver(SF::abstract_linear_solver< T, S > **sol)
Definition: SF_init.h:220
void compute_surface_mesh(const meshdata< T, S > &mesh, const SF_nbr numbering, const hashmap::unordered_set< T > &tags, meshdata< T, S > &surfmesh)
Compute the surface of a given mesh.
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:917
void make_global(const vector< T > &vec, vector< T > &out, MPI_Comm comm)
make a parallel vector global
Definition: SF_network.h:225
void unique_resize(vector< T > &_P)
Definition: SF_sort.h:348
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 local_petsc_to_nodal_mapping(const meshdata< T, S > &mesh, index_mapping< T > &petsc_to_nodal)
T local_nodal_to_local_petsc(const meshdata< T, S > &mesh, int rank, T local_nodal)
void assemble_lumped_matrix(abstract_matrix< T, S > &mat, meshdata< mesh_int_t, mesh_real_t > &domain, matrix_integrator< mesh_int_t, mesh_real_t > &integrator)
bool is_init(const abstract_vector< T, S > *v)
void layout_from_count(const T count, vector< T > &layout, MPI_Comm comm)
Definition: SF_network.h:201
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
@ 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
void dup_IMP_node_state(IonIfBase &IF, int from, int to, GlobalData_t **localdata)
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
void get_kappa(sf_vec &kappa, IMPregion *ir, limpet::MULTI_IF &miif, double k)
compute the vector
Definition: electrics.cc:809
physic_t
Identifier for the different physics we want to set up.
Definition: physics_types.h:51
int stimidx_from_timeridx(const SF::vector< stimulus > &stimuli, const int timer_id)
determine link between timer and stimulus
Definition: electrics.cc:781
@ iotm_chkpt_list
Definition: timer_utils.h:44
@ iotm_console
Definition: timer_utils.h:44
@ iotm_spacedt
Definition: timer_utils.h:44
@ iotm_trace
Definition: timer_utils.h:44
@ iotm_chkpt_intv
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:833
void sample_wave_form(stim_pulse &sp, int idx)
sample a signal given in analytic form
Definition: stimulate.cc:337
void read_el_scale_vec(const char *file, mesh_t mt, SF::vector< double > &el_scale, int &el_scale_dpn)
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
@ sum_cond
Definition: fem_types.h:43
@ intra_cond
Definition: fem_types.h:43
@ para_cond
Definition: fem_types.h:43
void print_act_log(FILE_SPEC logger, const SF::vector< Activation > &acts, int idx)
Definition: electrics.cc:1506
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.
bool is_dbc(stim_t type)
whether stimulus is a dirichlet type. implies boundary conditions on matrix
Definition: stimulate.cc:77
SF::meshdata< mesh_int_t, mesh_real_t > sf_mesh
Definition: sf_interface.h:48
@ constPulse
Definition: stimulate.h:74
void compute_restr_idx_async(sf_mesh &mesh, SF::vector< mesh_int_t > &inp_idx, SF::vector< mesh_int_t > &idx)
Definition: electrics.cc:577
void apply_stim_to_vector(const stimulus &s, sf_vec &vec, bool add)
Definition: electrics.cc:452
void recover_phie_std(sf_vec &vm, phie_recovery_data &rcv)
Definition: electrics.cc:1935
int set_dir(IO_t dest)
Definition: sim_utils.cc:495
@ ACT_THRESH
Definition: electrics.h:172
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:292
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
V dist(const vec3< V > &p1, const vec3< V > &p2)
Definition: vect.h:114
@ Phi_ex
Definition: stimulate.h:78
@ Vm_clmp
Definition: stimulate.h:78
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:401
void init_stim_info(void)
uses potential for stimulation
Definition: stimulate.cc:49
int output_all_activations(FILE_SPEC fp, int *ibuf, double *act_tbuf, int nlacts)
Definition: electrics.cc:1610
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:138
void savequit()
save state and quit simulator
Definition: sim_utils.cc:1703
bool have_dbc_stims(const SF::vector< stimulus > &stimuli)
return wheter any stimuli require dirichlet boundary conditions
Definition: electrics.cc:856
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
void setup_dataout(const int dataout, std::string dataout_vtx, mesh_t grid, SF::vector< mesh_int_t > *&restr, bool async)
Definition: electrics.cc:612
char * get_file_dir(const char *file)
Definition: sim_utils.cc:1332
@ POSTPROC
Definition: sim_utils.h:55
@ CURDIR
Definition: sim_utils.h:55
@ OUTPUT
Definition: sim_utils.h:55
void init_sv_gvec(gvec_data &GVs, limpet::MULTI_IF *miif, sf_vec &tmpl, igb_output_manager &output_manager)
Definition: ionics.cc:606
void assemble_sv_gvec(gvec_data &gvecs, limpet::MULTI_IF *miif)
Definition: ionics.cc:677
void constant_total_stimulus_current(SF::vector< stimulus > &stimuli, sf_mat &mass_i, sf_mat &mass_e, limpet::MULTI_IF *miif, FILE_SPEC logger)
Scales stimulus current to maintain constant total current across affected regions.
Definition: electrics.cc:2252
int postproc_recover_phie()
Definition: electrics.cc:2002
char * dupstr(const char *old_str)
Definition: basics.cc:44
void balance_electrode(elliptic_solver &ellip, SF::vector< stimulus > &stimuli, int balance_from, int balance_to)
Definition: electrics.cc:386
void set_elec_tissue_properties(MaterialType *mtype, Electrics::grid_t g, FILE_SPEC logger)
Fill the RegionSpec of an electrics grid with the associated inputs from the param parameters.
Definition: electrics.cc:108
void compute_restr_idx(sf_mesh &mesh, SF::vector< mesh_int_t > &inp_idx, SF::vector< mesh_int_t > &idx)
Definition: electrics.cc:544
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
@ extra_elec_msh
Definition: sf_interface.h:61
@ phie_recv_msh
Definition: sf_interface.h:70
@ intra_elec_msh
Definition: sf_interface.h:60
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
void setup_phie_recovery_data(phie_recovery_data &data)
Definition: electrics.cc:2085
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
Basic_physic * get_physics(physic_t p, bool error_if_missing)
Convinience function to get a physics.
Definition: sim_utils.cc:909
const char * get_tsav_ext(double time)
Definition: electrics.cc:865
SF::abstract_matrix< SF_int, SF_real > sf_mat
Definition: sf_interface.h:52
void compute_IIF(limpet::IonIfBase &pIF, limpet::GlobalData_t **impdata, int n)
Definition: ionics.cc:504
V timing(V &t2, const V &t1)
Definition: basics.h:464
void read_indices(SF::vector< T > &idx, const std::string filename, const hashmap::unordered_map< mesh_int_t, mesh_int_t > &dd_map, MPI_Comm comm)
Read indices from a file.
Definition: fem_utils.h:120
void update_cwd()
save the current working directory to curdir so that we can switch back to it if needed.
Definition: sim_utils.cc:490
void f_close(FILE_SPEC &f)
Close a FILE_SPEC.
Definition: basics.cc:165
@ ElecMat
Definition: fem_types.h:39
vec3< POINT_REAL > Point
Definition: vect.h:93
file_desc * FILE_SPEC
Definition: basics.h:138
#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 DATAOUT_SURF
Definition: sim_utils.h:60
#define BIDOMAIN
Definition: sim_utils.h:155
#define DATAOUT_VOL
Definition: sim_utils.h:61
#define MONODOMAIN
Definition: sim_utils.h:154
#define EXP_POSTPROCESS
Definition: sim_utils.h:172
#define DATAOUT_NONE
Definition: sim_utils.h:59
#define PSEUDO_BIDM
Definition: sim_utils.h:156
#define DATAOUT_VTX
Definition: sim_utils.h:62
Electrical stimulation functions.
SF_int niter
number of iterations
SF_int reason
number of iterations
std::string name
the solver name
virtual void setup_solver(abstract_matrix< T, S > &mat, double tol, int max_it, short norm, std::string name, bool has_nullspace, void *logger, const char *solver_opts_file, const char *default_opts)=0
event detection data structures
Definition: electrics.h:175
description of materal properties in a mesh
Definition: fem_types.h:121
SF::vector< RegionSpecs > regions
array with region params
Definition: fem_types.h:126
SF::vector< double > el_scale
optionally provided per-element params scale
Definition: fem_types.h:127
int el_scale_dpn
0=disabled, 1=isotropic scalar, 3=anisotropic (sl, st, sn) per element
Definition: fem_types.h:128
region based variations of arbitrary material parameters
Definition: fem_types.h:93
physMaterial * material
material parameter description
Definition: fem_types.h:98
int nsubregs
#subregions forming this region
Definition: fem_types.h:96
int * subregtags
FEM tags forming this region.
Definition: fem_types.h:97
char * regname
name of region
Definition: fem_types.h:94
int regID
region ID
Definition: fem_types.h:95
bool activated
flag sentinel activation
Definition: electrics.h:199
int ID
ID of LAT detector used as sentinel.
Definition: electrics.h:203
double t_start
start of observation window
Definition: electrics.h:200
double t_window
duration of observation window
Definition: electrics.h:201
double t_quiesc
measure current duration of quiescence
Definition: electrics.h:202
double ExVal[3]
extracellular conductivity eigenvalues
Definition: fem_types.h:62
cond_t g
rule to build conductivity tensor
Definition: fem_types.h:64
double InVal[3]
intracellular conductivity eigenvalues
Definition: fem_types.h:61
double BathVal[3]
bath conductivity eigenvalues
Definition: fem_types.h:63
File descriptor struct.
Definition: basics.h:133
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
void init_logger(const char *filename)
Definition: timers.cc:11
void log_stats(double tm, bool cflg)
Definition: timers.cc:27
void update_iter(const int curiter)
Definition: timers.cc:69
double slvtime
total solver time
Definition: timers.h:21
sf_vec * phie_rec
The phie recovery output vector buffer.
Definition: electrics.h:242
SF::vector< mesh_real_t > pts
The phie recovery locations.
Definition: electrics.h:241
physMat_t material_type
ID of physics material.
Definition: fem_types.h:53