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