openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
electrics_eikonal.cc
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #include "electrics_eikonal.h"
28 #include <cstdlib>
29 #include "SF_globals.h"
30 #include "basics.h"
31 #include "petsc_utils.h"
32 #include "timers.h"
33 #include "stimulate.h"
34 
35 #include "SF_init.h"
36 
37 namespace opencarp
38 {
39 
41 {
42  if (get_size() > 1) {
43  log_msg(NULL, 5, 0, "DREAM/Eikonal physics currently do not support MPI parallelization. Use openMP instead. Aborting!");
44  EXIT(EXIT_FAILURE);
45  }
46 
47  double t1, t2;
48  get_time(t1);
49 
50  set_dir(OUTPUT);
51 
52  // open logger
53  logger = f_open("eikonal.log", param_globals::experiment != 4 ? "w" : "r");
54 
55  // setup mappings between extra and intra grids, algebraic and nodal,
56  // and between PETSc and canonical orderings
57  setup_mappings();
58 
59  eik_tech = static_cast<Eikonal::eikonal_t>(param_globals::dream.solve);
60 
61  // the ionic physics is currently triggered from inside the Electrics to have tighter
62  // control over it. The standalone eikonal solver does not require it.
63  switch (eik_tech) {
64  case EIKONAL: break;
65  default:
66  ion.logger = logger;
67  ion.initialize();
68  }
69 
70  // set up Intracellular tissue
71  set_elec_tissue_properties(mtype, intra_grid, logger);
72  region_mask(intra_elec_msh, mtype[intra_grid].regions, mtype[intra_grid].regionIDs, true, "gregion_i");
73 
74  // add electrics timer for time stepping, add to time stepper tool (TS)
75  double global_time = user_globals::tm_manager->time;
76  timer_idx = user_globals::tm_manager->add_eq_timer(global_time, param_globals::tend, 0,
77  param_globals::dt, 0, "elec::ref_dt", "TS");
78 
79  // electrics stimuli setup
80  setup_stimuli();
81 
82  // set up the linear equation systems. this needs to happen after the stimuli have been
83  // set up, since we need boundary condition info
84  setup_solvers();
85 
86  // the next setup steps require the solvers to be set up, since they use the matrices
87  // generated by those
88 
89  // initialize the LATs detector
90  switch (eik_tech) {
91  case EIKONAL: break; // not available for pure eikonal solve
92  default:
94  }
95 
96  // prepare the electrics output. we skip it if we do post-processing
97  if (param_globals::experiment != EXP_POSTPROCESS)
98  setup_output();
99 
100  this->initialize_time += timing(t2, t1);
101  log_msg(NULL, 0, 0, "All done in %f sec.", float(t2 - t1));
102 }
103 
104 void Eikonal::set_elec_tissue_properties(MaterialType* mtype, Eikonal::grid_t g, FILE_SPEC logger)
105 {
106  MaterialType* m = mtype + g;
107 
108  // initialize random conductivity fluctuation structure with PrM values
109  m->regions.resize(param_globals::num_gregions);
110 
111  const char* grid_name = g == Eikonal::intra_grid ? "intracellular" : "extracellular";
112  log_msg(logger, 0, 0, "Setting up %s tissue poperties for %d regions ..", grid_name,
113  param_globals::num_gregions);
114 
115  char buf[64];
116  RegionSpecs* reg = m->regions.data();
117 
118  for (size_t i = 0; i < m->regions.size(); i++, reg++) {
119  if (!strcmp(param_globals::gregion[i].name, "")) {
120  snprintf(buf, sizeof buf, ", gregion_%d", int(i));
121  param_globals::gregion[i].name = dupstr(buf);
122  }
123 
124  reg->regname = strdup(param_globals::gregion[i].name);
125  reg->regID = i;
126  reg->nsubregs = param_globals::gregion[i].num_IDs;
127  if (!reg->nsubregs)
128  reg->subregtags = NULL;
129  else {
130  reg->subregtags = new int[reg->nsubregs];
131 
132  for (int j = 0; j < reg->nsubregs; j++)
133  reg->subregtags[j] = param_globals::gregion[i].ID[j];
134  }
135 
136  // describe material in given region
137  elecMaterial* emat = new elecMaterial();
138  emat->material_type = ElecMat;
139 
140  emat->InVal[0] = param_globals::gregion[i].g_il;
141  emat->InVal[1] = param_globals::gregion[i].g_it;
142  emat->InVal[2] = param_globals::gregion[i].g_in;
143 
144  emat->ExVal[0] = param_globals::gregion[i].g_el;
145  emat->ExVal[1] = param_globals::gregion[i].g_et;
146  emat->ExVal[2] = param_globals::gregion[i].g_en;
147 
148  emat->BathVal[0] = param_globals::gregion[i].g_bath;
149  emat->BathVal[1] = param_globals::gregion[i].g_bath;
150  emat->BathVal[2] = param_globals::gregion[i].g_bath;
151 
152  // convert units from S/m -> mS/um
153  for (int j = 0; j < 3; j++) {
154  emat->InVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
155  emat->ExVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
156  emat->BathVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
157  }
158  reg->material = emat;
159  }
160 
161  if ((g == Eikonal::intra_grid && strlen(param_globals::gi_scale_vec)) ||
162  (g == Eikonal::extra_grid && strlen(param_globals::ge_scale_vec))) {
164  sf_mesh& mesh = get_mesh(mt);
165  const char* file = g == Eikonal::intra_grid ? param_globals::gi_scale_vec : param_globals::ge_scale_vec;
166 
167  size_t num_file_entries = SF::root_count_ascii_lines(file, mesh.comm);
168 
169  if (num_file_entries != mesh.g_numelem)
170  log_msg(0, 4, 0, "%s warning: number of %s conductivity scaling entries does not match number of elements!",
171  __func__, get_mesh_type_name(mt));
172 
173  // set up parallel element vector and read data
174  sf_vec* escale;
175  SF::init_vector(&escale, get_mesh(mt), 1, sf_vec::elemwise);
176  escale->read_ascii(g == Eikonal::intra_grid ? param_globals::gi_scale_vec : param_globals::ge_scale_vec);
177 
178  if (get_size() > 1) {
179  // set up element vector permutation and permute
180  if (get_permutation(mt, ELEM_PETSC_TO_CANONICAL, 1) == NULL) {
182  }
184  sc(*escale, false);
185  }
186 
187  // copy data into SF::vector
188  SF_real* p = escale->ptr();
189  m->el_scale.assign(p, p + escale->lsize());
190  escale->release_ptr(p);
191  }
192 }
193 
194 void Eikonal::setup_mappings()
195 {
196  bool intra_exits = mesh_is_registered(intra_elec_msh), extra_exists = mesh_is_registered(extra_elec_msh);
197  assert(intra_exits);
198  const int dpn = 1;
199 
200  // It may be that another physic (e.g. ionic models) has already computed the intracellular mappings,
201  // thus we first test their existence
202  if (get_scattering(intra_elec_msh, ALG_TO_NODAL, dpn) == NULL) {
203  log_msg(logger, 0, 0, "%s: Setting up intracellular algebraic-to-nodal scattering.", __func__);
205  }
207  log_msg(logger, 0, 0, "%s: Setting up intracellular PETSc to canonical permutation.", __func__);
209  }
210 }
211 
213 {
214  double t1, t2;
215  get_time(t1);
216 
217  switch (eik_tech) {
218  case EIKONAL: solve_EIKONAL(); break;
219  case DREAM: solve_DREAM(); break;
220  default: solve_RE(); break;
221  }
222 
224  // output lin solver stats
226  }
227  this->compute_time += timing(t2, t1);
228 
229  // since the traces have their own timing, we check for trace dumps in the compute step loop
232 }
233 
235 {
236  double t1, t2;
237  get_time(t1);
238 
239  switch (eik_tech) {
240  case EIKONAL: break;
241  default:
243  output_manager_time.write_data(); // does not exist in EIKONAL
244  }
245 
246  if (do_output_eikonal) {
248  do_output_eikonal = false;
249  }
250 
251  double curtime = timing(t2, t1);
252  this->output_time += curtime;
253 
254  IO_stats.calls++;
255  IO_stats.tot_time += curtime;
256 
259 }
260 
265 {
266  // close logger
267  f_close(logger);
268 
269  switch (eik_tech) {
270  case EIKONAL:
272  break;
273  default:
274  // output LAT data
278  ion.destroy();
279  }
280 }
281 
282 void Eikonal::setup_stimuli()
283 {
284  // initialize basic stim info data (used units, supported types, etc)
285  init_stim_info();
286 
287  stimuli.resize(param_globals::num_stim);
288  for (int i = 0; i < param_globals::num_stim; i++) {
289  // construct new stimulus
290  stimulus& s = stimuli[i];
291 
292  if (param_globals::stim[i].crct.type != 0 && param_globals::stim[i].crct.type != 9) {
293  // In the Eikonal class, only the intracellular domain is registered.
294  // Therefore, only I_tm and Vm_clmp are compatible.
295  log_msg(NULL, 5, 0, "%s error: stimulus of type %i is incompatible with the eikonal model! Use I_tm or Vm_clmp instead. Aborting!", __func__, s.phys.type);
296  EXIT(EXIT_FAILURE);
297  }
298 
300  s.translate(i);
301 
302  s.setup(i);
303 
304  if (s.electrode.dump_vtx)
305  s.dump_vtx_file(i);
306 
307  log_msg(NULL, 2, 0, "Only geometry, start time, npls, and bcl of stim[%i] are used", i);
308 
309  if (param_globals::stim[i].pulse.dumpTrace && get_rank() == 0) {
310  set_dir(OUTPUT);
311  s.pulse.wave.write_trace(s.name + ".trc");
312  }
313  }
314 
316 }
317 
318 void Eikonal::stimulate_intracellular()
319 {
320  parabolic_solver& ps = parab_solver;
321 
322  // iterate over stimuli
323  for (stimulus& s : stimuli) {
324  if (s.is_active()) {
325  // for active stimuli, deal with the stimuli-type specific stimulus application
326  switch (s.phys.type) {
327  case I_tm: {
328  if (param_globals::operator_splitting) {
329  apply_stim_to_vector(s, *ps.Vmv, true);
330  } else {
331  SF_real Cm = 1.0;
332  timer_manager& tm = *user_globals::tm_manager;
333  SF_real sc = tm.time_step / Cm;
334 
335  ps.Irhs->set(0.0);
336  apply_stim_to_vector(s, *ps.Irhs, true);
337 
338  *ps.tmp_i1 = *ps.IIon;
339  *ps.tmp_i1 -= *ps.Irhs;
340  *ps.tmp_i1 *= sc; // tmp_i1 = sc * (IIon - Irhs)
341 
342  // add ionic, transmembrane and intracellular currents to rhs
343  if (param_globals::parab_solve != parabolic_solver::EXPLICIT)
344  ps.mass_i->mult(*ps.tmp_i1, *ps.Irhs);
345  else
346  *ps.Irhs = *ps.tmp_i1;
347  }
348  break;
349  }
350 
351  case Illum: {
352  sf_vec* illum_vec = ion.miif->gdata[limpet::illum];
353 
354  if (illum_vec == NULL) {
355  log_msg(0, 5, 0, "Cannot apply illumination stim: global vector not present!");
356  EXIT(EXIT_FAILURE);
357  } else {
358  apply_stim_to_vector(s, *illum_vec, false);
359  }
360 
361  break;
362  }
363 
364  default: break;
365  }
366  }
367  }
368 }
369 
370 void Eikonal::clamp_Vm()
371 {
372  for (stimulus& s : stimuli) {
373  if (s.phys.type == Vm_clmp && s.is_active())
375  }
376 }
377 
378 void Eikonal::setup_output()
379 {
380  int rank = get_rank();
381  SF::vector<mesh_int_t>* restr_i = NULL;
382  SF::vector<mesh_int_t>* restr_e = NULL;
383  set_dir(OUTPUT);
384 
385  setup_dataout(param_globals::dataout_i, param_globals::dataout_i_vtx, intra_elec_msh,
386  restr_i, param_globals::num_io_nodes > 0);
387 
388  if (param_globals::dataout_i) {
389  switch (eik_tech) {
390  case EIKONAL:
391  output_manager_cycle.register_output(eik_solver.AT, intra_elec_msh, 1, param_globals::dream.output.atfile, "ms", restr_i);
392  break;
393  case DREAM:
394  output_manager_time.register_output(parab_solver.Vmv, intra_elec_msh, 1, param_globals::vofile, "mV", restr_i);
395  output_manager_cycle.register_output(eik_solver.AT, intra_elec_msh, 1, param_globals::dream.output.atfile, "ms", restr_i);
396  output_manager_cycle.register_output(eik_solver.RT, intra_elec_msh, 1, param_globals::dream.output.rtfile, "ms", restr_i);
397  if (strcmp(param_globals::dream.output.idifffile, "") != 0) {
398  output_manager_time.register_output(eik_solver.Idiff, intra_elec_msh, 1, param_globals::dream.output.idifffile, "muA/cm2", restr_i);
399  }
400  break;
401  default:
402  output_manager_time.register_output(parab_solver.Vmv, intra_elec_msh, 1, param_globals::vofile, "mV", restr_i);
403  output_manager_cycle.register_output(eik_solver.AT, intra_elec_msh, 1, param_globals::dream.output.atfile, "ms", restr_i);
404  if (strcmp(param_globals::dream.output.idifffile, "") != 0) {
405  output_manager_time.register_output(eik_solver.Idiff, intra_elec_msh, 1, param_globals::dream.output.idifffile, "muA/cm2", restr_i);
406  }
407  }
408  }
409 
411 
412  if (param_globals::num_trace) {
413  sf_mesh& imesh = get_mesh(intra_elec_msh);
414  open_trace(ion.miif, param_globals::num_trace, param_globals::trace_node, NULL, &imesh);
415  }
416 
417  // initialize generic logger for IO timings per time_dt
418  IO_stats.init_logger("IO_stats.dat");
419 }
420 
421 void Eikonal::dump_matrices()
422 {
423  std::string bsname = param_globals::dump_basename;
424  std::string fn;
425 
426  set_dir(OUTPUT);
427 
428  // dump monodomain matrices
429  if (param_globals::parab_solve == 1) {
430  // using Crank-Nicolson
431  fn = bsname + "_Ki_CN.bin";
432  parab_solver.lhs_parab->write(fn.c_str());
433  }
434  fn = bsname + "_Ki.bin";
435  parab_solver.rhs_parab->write(fn.c_str());
436 
437  fn = bsname + "_Mi.bin";
438  parab_solver.mass_i->write(fn.c_str());
439 }
440 
443 double Eikonal::timer_val(const int timer_id)
444 {
445  // determine
446  int sidx = stimidx_from_timeridx(stimuli, timer_id);
447  double val = 0.0;
448  if (sidx != -1) {
449  stimuli[sidx].value(val);
450  } else
451  val = std::nan("NaN");
452 
453  return val;
454 }
455 
458 std::string Eikonal::timer_unit(const int timer_id)
459 {
460  int sidx = stimidx_from_timeridx(stimuli, timer_id);
461  std::string s_unit;
462 
463  if (sidx != -1)
464  // found a timer-linked stimulus
465  s_unit = stimuli[sidx].pulse.wave.f_unit;
466 
467  return s_unit;
468 }
469 
470 void Eikonal::setup_solvers()
471 {
472  set_dir(OUTPUT);
473 
474  switch (eik_tech) {
475  case EIKONAL:
476  eik_solver.init();
477  break;
478  default:
479  parab_solver.init();
481  eik_solver.init();
483  if (param_globals::dump2MatLab) dump_matrices();
484  }
485 }
486 
487 void Eikonal::checkpointing()
488 {
489  const timer_manager& tm = *user_globals::tm_manager;
490 
491  // regular user selected state save
492  if (tm.trigger(iotm_chkpt_list)) {
493  char save_fnm[1024];
494  const char* tsav_ext = get_tsav_ext(tm.time);
495 
496  snprintf(save_fnm, sizeof save_fnm, "%s.%s.roe", param_globals::write_statef, tsav_ext);
497 
498  ion.miif->dump_state(save_fnm, tm.time, intra_elec_msh, false, GIT_COMMIT_COUNT);
499  eik_solver.save_eikonal_state(tsav_ext);
500  }
501 
502  // checkpointing based on interval
503  if (tm.trigger(iotm_chkpt_intv)) {
504  char save_fnm[1024];
505  snprintf(save_fnm, sizeof save_fnm, "checkpoint.%.1f.roe", tm.time);
506  ion.miif->dump_state(save_fnm, tm.time, intra_elec_msh, false, GIT_COMMIT_COUNT);
507  }
508 }
509 
510 void Eikonal::solve_EIKONAL()
511 {
512  if (user_globals::tm_manager->time == 0) {
513  eik_solver.FIM();
516  do_output_eikonal = true;
517  }
518 }
519 
520 void Eikonal::solve_RE()
521 {
522  if (user_globals::tm_manager->time == 0) {
523  eik_solver.FIM();
526  do_output_eikonal = true;
527  }
528  solve_RD();
529 }
530 
531 void Eikonal::solve_DREAM()
532 {
534  solve_RD();
535  } else {
536  if (user_globals::tm_manager->time != 0) {
538  }
539 
540  eik_solver.cycFIM();
542 
545  do_output_eikonal = true;
546  }
547 }
548 
549 void Eikonal::solve_RD()
550 {
551  // if requested, we checkpoint the current state
552  checkpointing();
553 
554  // activation checking
555  const double time = user_globals::tm_manager->time,
556  time_step = user_globals::tm_manager->time_step;
557 
558  lat.check_acts(time);
559  lat.check_quiescence(time, time_step);
560 
561  clamp_Vm();
562 
563  *parab_solver.old_vm = *parab_solver.Vmv; // needed in step D of DREAM
564 
565  // compute ionics update
566  ion.compute_step();
567 
568  // Compute the I diff current
569  *eik_solver.Idiff *= 0;
572 
573  if (eik_tech == REp) {
575  }
576 
577  switch (eik_tech) {
579  default: break;
580  }
581 
582  clamp_Vm();
583 }
584 
586 {
587  if (param_globals::output_level > 1) log_msg(0, 0, 0, "\n *** Initializing Eikonal Solver ***\n");
588  stats.init_logger("eik_stats.dat");
589 
590  if (param_globals::dream.output.debugNode >= 0) {
591  nodeData.idX = param_globals::dream.output.debugNode;
592  char buf[256];
593  snprintf(buf, sizeof buf, "node_%d.dat", nodeData.idX);
594  nodeData.init_logger(buf);
595  }
596 
597  // currently only used for output
598  const sf_mesh& mesh = get_mesh(intra_elec_msh);
599  twoFib = mesh.she.size() > 0;
603 
604  num_pts = mesh.l_numpts;
605 
606  e2n_con = mesh.con; // Connectivity vector with nodes that belong to each element
607  // The 4 nodes from index 4j to 4j+3 belong to the same element for 0<=j<num elements
608  elem_start = mesh.dsp; // For the i_th element elem_start[j] stores its starting position in the "connect" vector
609 
613 
614  List.assign(num_pts, 0); // List of active nodes
615  T_A.assign(num_pts, inf); // Activation times
616  D_I.assign(num_pts, 0); // Diastolic Intervals
617  T_R.assign(num_pts, -400); // Recovery times
618  TA_old.assign(num_pts, inf); // Activation Time is the previous cycle
619  num_changes.assign(num_pts, 0); // Number of Times entered in the List per current LAT
620  nReadded2List.assign(num_pts, 0); // Number of Times entered in the List per current LAT
621  stim_status.assign(num_pts, 0); // Status of nodes
622 
623  switch (mesh.type[0]) {
624  case 0:
625  MESH_SIZE = 4;
626  break;
627 
628  case 6:
629  MESH_SIZE = 3;
630  break;
631 
632  case 7:
633  MESH_SIZE = 2;
634  break;
635 
636  default:
637  log_msg(0, 5, 0, "Error: Type of element is not compatible with this version of the eikonal model. Use tetrahedra or triangles");
638  EXIT(EXIT_FAILURE);
639  break;
640  }
641 
642  if (strlen(param_globals::start_statef) > 0) load_state_file();
643 
644  create_node_to_node_graph();
645 
646  translate_stim_to_eikonal();
647 
648  precompute_squared_anisotropy_metric();
649 }
650 
652 {
653  double t1, t2;
654  get_time(t1);
655 
656  const sf_mesh& mesh = get_mesh(intra_elec_msh);
657 
658  diff_cur.resize(mesh.l_numpts);
659  rho_cvrest.resize(mesh.l_numpts);
663 
664  // --- Precompute resolved region IDs for each point ---
665  // Use the same region delegation as the ionics class
666  SF::vector<int> regionIDs;
667  SF::vector<RegionSpecs> rs(param_globals::num_imp_regions);
668  for (size_t i = 0; i < rs.size(); i++) {
669  rs[i].nsubregs = param_globals::imp_region[i].num_IDs;
670  rs[i].subregtags = param_globals::imp_region[i].ID;
671  for (int j = 0; j < rs[i].nsubregs; j++) {
672  if (rs[i].subregtags[j] == -1 && get_rank() == 0)
673  log_msg(NULL, 3, ECHO, "Warning: not all %u IDs provided for imp_region[%u]!\n", rs[i].nsubregs, i);
674  }
675  }
676  if (rs.size() == 1) {
677  regionIDs.assign(mesh.l_numpts, 0);
678  } else {
679  region_mask(intra_elec_msh, rs, regionIDs, false, "imp_regions");
680  }
681 
682  // --- Parallel initialization, one write per vertex ---
683  #pragma omp parallel for schedule(dynamic)
684  for (int v = 0; v < mesh.l_numpts; v++) {
685  int reg = regionIDs[v];
686 
687  const auto& region_diff = param_globals::imp_region[reg].dream.Idiff;
688  const auto& region_rest = param_globals::imp_region[reg].dream.CVrest;
689 
690  auto model = static_cast<eikonal_solver::Idiff_t>(region_diff.model);
691  diff_cur[v].model = model;
692 
693  if (model == GAUSS) {
694  diff_cur[v].alpha_1 = region_diff.alpha_i[0];
695  diff_cur[v].alpha_2 = region_diff.alpha_i[1];
696  diff_cur[v].alpha_3 = region_diff.alpha_i[2];
697  diff_cur[v].beta_1 = region_diff.beta_i[0];
698  diff_cur[v].beta_2 = region_diff.beta_i[1];
699  diff_cur[v].beta_3 = region_diff.beta_i[2];
700  diff_cur[v].gamma_1 = region_diff.gamma_i[0];
701  diff_cur[v].gamma_2 = region_diff.gamma_i[1];
702  diff_cur[v].gamma_3 = region_diff.gamma_i[2];
703  } else {
704  diff_cur[v].A_F = region_diff.A_F;
705  diff_cur[v].tau_F = region_diff.tau_F;
706  diff_cur[v].V_th = region_diff.V_th;
707  }
708 
709  rho_cvrest[v] = region_rest.rho;
710  kappa_cvrest[v] = region_rest.kappa;
711  theta_cvrest[v] = region_rest.theta;
712  denom_cvrest[v] = log(region_rest.rho) / region_rest.psi;
713  }
714  if (param_globals::output_level) log_msg(NULL, 0, 0, "Diffusion current and CV restitution initialized in %f sec.", timing(t2, t1));
715 }
716 
717 void eikonal_solver::translate_stim_to_eikonal()
718 {
719  double t1, t2;
720  get_time(t1);
721  Index_currStim = 0;
722 
723  // Collect all pulses as (start_time, stimulus_index)
724  std::vector<std::pair<SF_real, int>> stim_events;
725  stim_events.reserve(param_globals::num_stim * 8); // rough guess
726 
727  for (int stim_idx = 0; stim_idx < stimuliRef->size(); ++stim_idx) {
728  const stimulus& s = (*stimuliRef)[stim_idx];
729  for (int idx_pls = 0; idx_pls < s.ptcl.npls; ++idx_pls) {
730  SF_real start_time = s.ptcl.start + s.ptcl.pcl * idx_pls;
731  stim_events.emplace_back(start_time, stim_idx);
732  }
733  }
734 
735  // Sort events by start_time
736  std::sort(stim_events.begin(), stim_events.end(),
737  [](auto& a, auto& b) { return a.first < b.first; });
738 
739  // Reserve enough space for output
740  size_t total_nodes = 0;
741  for (auto& ev : stim_events)
742  total_nodes += (*stimuliRef)[ev.second].electrode.vertices.size();
743  StimulusPoints.resize(total_nodes + 1, -1);
744  StimulusTimes.resize(total_nodes + 1, -1);
745 
746  // Fill outputs
747  size_t count = 0;
748  for (auto& ev : stim_events) {
749  const stimulus& s = (*stimuliRef)[ev.second];
750  for (mesh_int_t v : s.electrode.vertices) {
751  StimulusPoints[count] = v;
752  StimulusTimes[count] = ev.first;
753  ++count;
754  }
755  }
756 
757  if (param_globals::output_level)
758  log_msg(NULL, 0, 0, "Translating stimuli for eikonal model done in %f sec.", timing(t2, t1));
759 }
760 
761 void eikonal_solver::create_node_to_node_graph()
762 {
763  double t1, t2;
764  get_time(t1);
765 
766  n2n_dsp.resize(num_pts + 1, 0);
767  const sf_mesh& mesh = get_mesh(intra_elec_msh);
768 
769  // Thread-local storage for neighbors
770  std::vector<std::vector<mesh_int_t>> thread_neighbors(num_pts);
771 
772  // Preallocate thread-local marker arrays for O(n) duplicate removal
773  #pragma omp parallel
774  {
775  std::vector<char> mark(num_pts, 0); // one per thread
776 
777  #pragma omp for schedule(dynamic)
778  for (int point_idx = 0; point_idx < num_pts; point_idx++) {
779  int numNBElem = n2e_dsp[point_idx + 1] - n2e_dsp[point_idx];
780  std::vector<mesh_int_t> neighbors;
781  neighbors.reserve(numNBElem * MESH_SIZE); // rough estimate
782 
783  // Loop over all elements containing this node
784  for (int eedsp = 0; eedsp < numNBElem; eedsp++) {
785  int currElem = n2e_con[n2e_dsp[point_idx] + eedsp];
786 
787  // Loop over all nodes of the current element
788  for (int nndsp = 0; nndsp < MESH_SIZE; nndsp++) {
789  int nb = e2n_con[elem_start[currElem] + nndsp];
790 
791  if (nb == point_idx) continue; // skip self
792  if (!mark[nb]) {
793  mark[nb] = 1;
794  neighbors.push_back(nb);
795  }
796  }
797  }
798 
799  // Reset marks for the next iteration
800  for (int nb : neighbors) mark[nb] = 0;
801 
802  // Move neighbor list into the final container
803  thread_neighbors[point_idx] = std::move(neighbors);
804  }
805  }
806 
807  // Build prefix sums (serial, cheap compared to above)
808  for (int i = 0; i < num_pts; i++) {
809  n2n_dsp[i + 1] = n2n_dsp[i] + thread_neighbors[i].size();
810  }
811 
812  // Allocate once, then fill in parallel
813  n2n_connect.resize(n2n_dsp[num_pts]);
814 
815  #pragma omp parallel for schedule(static)
816  for (int i = 0; i < num_pts; i++) {
817  std::copy(thread_neighbors[i].begin(),
818  thread_neighbors[i].end(),
819  n2n_connect.begin() + n2n_dsp[i]);
820  }
821 
822  if (param_globals::output_level) {
823  log_msg(NULL, 0, 0, "Node-to-node graph done in %f sec.", timing(t2, t1));
824  }
825 }
826 
827 void eikonal_solver::precompute_squared_anisotropy_metric()
828 {
829  double t1, t2;
830  get_time(t1);
831 
832  const sf_mesh& mesh = get_mesh(intra_elec_msh);
833 
834  S.resize(mesh.l_numelem); // store anisotropy matrices (dimensionless)
835  CV_L.resize(mesh.l_numelem); // store CV_L per element for later use
836 
837  // temporary variables
838  SF::dmat<double> I(3, 3);
839  I.assign(0.0);
840  I.diag(1.0);
841 
842  #pragma omp parallel
843  {
844  SF::dmat<double> Saniso(3, 3);
845  SF::dmat<double> Ff(3, 3);
846  SF::dmat<double> diff(3, 3);
847  SF::Point f, s, n;
848 
849  #pragma omp for schedule(static)
850  for (int eidx = 0; eidx < mesh.l_numelem; eidx++) {
851  double vl = 0.0, vt = 0.0, vn = 0.0;
852 
853  // Find region match (break early when found)
854  for (int g = 0; g < param_globals::num_gregions; g++) {
855  const auto& reg = param_globals::gregion[g];
856  const auto& dream = reg.dream;
857  for (int j = 0; j < reg.num_IDs; j++) {
858  if (mesh.tag[eidx] == reg.ID[j]) {
859  vl = dream.vel_l;
860  vt = dream.vel_t;
861  vn = dream.vel_n;
862  goto region_found; // break out of both loops
863  }
864  }
865  }
866  region_found:;
867 
868  // Store CV_L directly (for later rescaling)
869  CV_L[eidx] = vl;
870 
871  // Compute anisotropy ratios relative to CV_L
872  double AR_T2 = (vl / vt) * (vl / vt);
873  double AR_N2 = (vl / vn) * (vl / vn);
874 
875  f.x = mesh.fib[3 * eidx + 0];
876  f.y = mesh.fib[3 * eidx + 1];
877  f.z = mesh.fib[3 * eidx + 2];
878 
879  if (twoFib) {
880  s.x = mesh.she[3 * eidx + 0];
881  s.y = mesh.she[3 * eidx + 1];
882  s.z = mesh.she[3 * eidx + 2];
883  n = cross(f, s);
884 
885  SF::outer_prod(f, f, 1.0, Saniso[0], false);
886  SF::outer_prod(s, s, AR_T2, Saniso[0], true);
887  SF::outer_prod(n, n, AR_N2, Saniso[0], true);
888  } else {
889  SF::outer_prod(f, f, 1.0, Ff[0], false);
890  diff = I - Ff;
891  diff *= AR_T2;
892 
893  Saniso = Ff + diff;
894  }
895 
896  S[eidx] = Saniso;
897  }
898  }
899 
900  if (param_globals::output_level)
901  log_msg(NULL, 0, 0, "Anisotropy tensors precomputed in %f sec.", timing(t2, t1));
902 }
903 
905 {
906  double t0, t1;
907  get_time(t0);
908 
909  // --- 0) Reset all activation times to infinity
910  std::fill(T_A.begin(), T_A.end(), inf);
911 
912  // --- 1) Build initial active list = neighbors of stimuli
913  std::vector<mesh_int_t> activeList;
914  activeList.reserve(num_pts / 10); // heuristic reserve to avoid frequent reallocs
915  std::vector<char> in_active(num_pts, 0);
916 
917  // Mark stimuli
918  for (size_t si = 0; si < StimulusTimes.size() - 1; ++si) { // skip last entry, since its a -1 placeholder currently still used in DREAM
919  mesh_int_t s = StimulusPoints[si];
920  T_A[s] = StimulusTimes[si];
921  }
922  // Seed neighbors; separate loop to avoid adding stim nodes into the list if they are a neighbor
923  for (size_t si = 0; si < StimulusTimes.size() - 1; ++si) { // skip last entry, since its a -1 placeholder currently still used in DREAM
924  mesh_int_t s = StimulusPoints[si];
925  for (int off = n2n_dsp[s]; off < n2n_dsp[s + 1]; ++off) {
926  mesh_int_t nb = n2n_connect[off];
927  if (T_A[nb] == inf) {
928  add_to_active(activeList, in_active, nb);
929  }
930  }
931  }
932 
933  // --- 2) Iterative solve of activeList
934  int niter = 0;
935  while (!activeList.empty() && niter <= param_globals::dream.fim.max_iter) {
936  const std::vector<mesh_int_t> activeVec = std::move(activeList);
937  activeList.clear();
938 
939  #pragma omp parallel
940  {
941  std::vector<mesh_int_t> local_active;
942  local_active.reserve(64); // small local buffer
943 
944  #pragma omp for schedule(dynamic)
945  for (size_t i = 0; i < activeVec.size(); ++i) {
946  mesh_int_t id = activeVec[i];
947  remove_from_active(in_active, id);
948 
949  SF_real p = T_A[id];
950  SF_real q = update(id);
951 
952  #pragma omp atomic write
953  T_A[id] = q;
954 
955  if (std::fabs(p - q) < param_globals::dream.fim.tol) {
956  for (int off = n2n_dsp[id]; off < n2n_dsp[id + 1]; ++off) {
957  mesh_int_t nb = n2n_connect[off];
958  p = T_A[nb];
959  q = update(nb);
960  if (p > q) {
961  #pragma omp atomic write
962  T_A[nb] = q;
963  local_active.push_back(nb);
964  }
965  }
966  } else {
967  local_active.push_back(id); // non-converged -> add back
968  }
969  }
970 
971  // Merge local results
972  #pragma omp critical
973  {
974  for (mesh_int_t nb : local_active) {
975  add_to_active(activeList, in_active, nb);
976  }
977  }
978  }
979  ++niter;
980  }
981 
982  // collect iteration stats
983  stats.update_iter(niter);
984  auto [minIt, maxIt] = std::minmax_element(T_A.begin(), T_A.end());
985  actMIN = *minIt;
986  actMAX = *maxIt;
987 
988  // --- 3) Copy into AT for output
989  double* atc = AT->ptr();
990  const SF_real* t_a = T_A.data();
991  #pragma omp parallel for simd
992  for (mesh_int_t i = 0; i < num_pts; ++i)
993  atc[i] = (t_a[i] == inf ? -1.0 : t_a[i]);
994  AT->release_ptr(atc);
995 
996  // --- 4) Timing & list‐size stats
997  auto dur = timing(t1, t0);
998  stats.slvtime_A += dur;
999  stats.minAT = actMIN;
1000  stats.maxAT = actMAX;
1001  stats.activeList = activeList.size();
1002  stats.bc_status = true;
1004 }
1005 
1007 {
1008  int niter = 0;
1009  double t1, t0;
1010  get_time(t0);
1011 
1012  if (sum(List) == 0) {
1013  compute_bc();
1014  }
1015 
1016  double time2stop_eikonal = param_globals::dream.tau_inc;
1017  float maxadvance = user_globals::tm_manager->time + param_globals::dream.tau_s + param_globals::dream.tau_inc + param_globals::dream.tau_max;
1018 
1019  do {
1020  if (StimulusTimes[Index_currStim] <= maxadvance) {
1021  compute_bc();
1022  }
1023 
1024  for (mesh_int_t indX = 0; indX < List.size(); indX++) {
1025  SF_real p = T_A[indX];
1026  SF_real q;
1027 
1028  if (List[indX] == 0) continue;
1029 
1030  q = compute_coherence(indX);
1031 
1032  T_A[indX] = q;
1033  num_changes[indX] = num_changes[indX] + 1;
1034 
1035  if (q > maxadvance) continue;
1036 
1037  if (abs(p - q) < param_globals::dream.fim.tol || (num_changes[indX] > param_globals::dream.fim.max_iter) || (q == inf || p == inf)) {
1038  for (int ii = n2n_dsp[indX]; ii < n2n_dsp[indX + 1]; ii++) {
1039  mesh_int_t indXNB = n2n_connect[ii];
1040 
1041  if (List[indXNB] == 1) {
1042  continue;
1043  }
1044 
1045  SF_real pNB = T_A[indXNB];
1046  SF_real qNB;
1047 
1048  qNB = compute_coherence(indXNB);
1049 
1050  bool node_is_valid = add_node_neighbor_to_list(T_R[indXNB], pNB, qNB) && qNB != inf && qNB > user_globals::tm_manager->time;
1051  // This second condition is there to be able to add a node to the list when a new valid activation time
1052  // is found but a reentry would be blocked by the L2 parameter. Normally this condition does not make or brake the
1053  // simulation but would leave individual nodes inactivated, which is not ideal.
1054  bool ignore_L2_if_valid = node_is_valid && qNB > pNB && nReadded2List[indXNB] >= param_globals::dream.fim.max_addpt;
1055 
1056  if (node_is_valid && (nReadded2List[indXNB] < param_globals::dream.fim.max_addpt) || ignore_L2_if_valid) {
1057  if (qNB > pNB) {
1058  nReadded2List[indXNB] = 0;
1059  }
1060  T_A[indXNB] = qNB;
1061  nReadded2List[indXNB]++;
1062  num_changes[indXNB] = 0;
1063  List[indXNB] = 1;
1064  if (param_globals::dream.output.debugNode == indXNB) {
1066  nodeData.idXNB = indX;
1067  nodeData.nbn_T_A = q;
1068  }
1069  }
1070  }
1071 
1072  List[indX] = 0;
1073  if (param_globals::dream.output.debugNode == indX) {
1075  }
1076  }
1077  }
1078 
1079  SF_real actMIN_old = actMIN;
1080  SF_real progress_time;
1081 
1082  update_Ta_in_active_list();
1083 
1084  if (actMIN > actMIN_old) {
1085  progress_time = actMIN - actMIN_old;
1086  } else {
1087  progress_time = 0;
1088  }
1089 
1090  time2stop_eikonal -= progress_time;
1091  niter++;
1092 
1093  } while (time2stop_eikonal > 0 && sum(List) > 0);
1094 
1095  if (sum(List) == 0) {
1097  }
1098 
1099  // copy for igb output
1100  double* atc = AT->ptr();
1101  double* rpt = RT->ptr();
1102  for (mesh_int_t i = 0; i < List.size(); i++) {
1103  if (T_A[i] == inf) {
1104  // for better visualization in meshalyzer
1105  atc[i] = -1;
1106  } else {
1107  atc[i] = T_A[i];
1108  }
1109  rpt[i] = T_R[i];
1110  }
1111  AT->release_ptr(atc);
1112  RT->release_ptr(rpt);
1113 
1114  // treat solver statistics
1115  auto dur = timing(t1, t0);
1116  stats.slvtime_A += dur;
1117  stats.update_iter(niter);
1118  stats.minAT = actMIN;
1119  stats.maxAT = actMAX;
1120  stats.activeList = sum(List);
1122 
1123  // treat node stats
1124  if (param_globals::dream.output.debugNode >= 0) {
1128  }
1129 
1130 } // close iterate list
1131 
1132 SF_real eikonal_solver::update(mesh_int_t& indX, SF_real CVrest_factor, bool isDREAM)
1133 {
1134  switch (MESH_SIZE) {
1135  case 4: return update_impl<4>(indX, CVrest_factor, isDREAM);
1136  case 3: return update_impl<3>(indX, CVrest_factor, isDREAM);
1137  case 2: return update_impl<2>(indX, CVrest_factor, isDREAM);
1138  default:
1139  return T_A[indX]; // fallback
1140  }
1141 }
1142 
1143 template <int N>
1144 SF_real eikonal_solver::update_impl(mesh_int_t& indX, SF_real CVrest_factor, bool CheckValidity)
1145 {
1146  double min = inf;
1147  const double time = user_globals::tm_manager->time;
1148  const sf_mesh& mesh = get_mesh(intra_elec_msh);
1149 
1150  // element-wise update
1151  for (int e_i = n2e_dsp[indX]; e_i < n2e_dsp[indX + 1]; e_i++) { // gather all elements the node belongs to
1152  int Elem_i = n2e_con[e_i];
1153  mesh_int_t indEle = elem_start[Elem_i];
1154  std::array<SF::Point, N> base; // points
1155  std::array<double, N> values; // activation times
1156  std::array<int, N> nodeIDs;
1157 
1158  std::size_t k = 0;
1159  for (std::size_t j = 0; j < N; j++) { // loop over nodes of element e_i
1160  int n_i = e2n_con[indEle + j];
1161  if (n_i != indX) {
1162  base[k].x = mesh.xyz[3 * n_i + 0]; base[k].y = mesh.xyz[3 * n_i + 1]; base[k].z = mesh.xyz[3 * n_i + 2];
1163  values[k] = T_A[n_i];
1164  nodeIDs[k] = n_i;
1165  k++;
1166  } else { // last slot is reserved for the current vertex we are solving for
1167  base[N - 1].x = mesh.xyz[3 * indX + 0]; base[N - 1].y = mesh.xyz[3 * indX + 1]; base[N - 1].z = mesh.xyz[3 * indX + 2];
1168  values[N - 1] = T_A[indX];
1169  nodeIDs[N - 1] = indX;
1170  }
1171  }
1172  // scale slowness metric by CV
1173  double cv = CV_L[Elem_i] * CVrest_factor;
1174  if (cv == 0.0) continue; // avoid 1/(cv*cv)
1175  SF::dmat<double> D = 1.0 / (cv * cv) * S[Elem_i];
1176 
1177  // 1) solve full N-simplex (eikonal run OR if valid for DREAM)
1178  if (!CheckValidity || !is_not_valid_update<N>(nodeIDs, time)) {
1179  LocalSolver<N> solver(D, base, values);
1180  SF_real tmp = solver.solve();
1181  if (min > tmp && tmp > T_R[indX] && compute_H(indX, tmp) > 0.0) {
1182  min = tmp;
1183  }
1184  continue;
1185  }
1186 
1187  // If the full N-simplex was not valid, we have to do additional checks for the DREAM,
1188  // since a subsimplex could be valid if e.g. only one node of a tet/triangle is invalid
1189  // 2) triangles that include the target (only done if MESH_SIZE is 4)
1190  if constexpr (N - 1 == 3) {
1191  // neighbor indices are 0..(N-2); choose pairs (i,j) and add target (N-1)
1192  for (int i = 0; i < (N - 1); ++i) {
1193  for (int j = i + 1; j < (N - 1); ++j) {
1194  // build nodeIDs/points/values for face {i, j, target}
1195  std::array<int, 3> tri_ids{nodeIDs[i], nodeIDs[j], nodeIDs[N - 1]};
1196  if (is_not_valid_update<3>(tri_ids, time)) continue;
1197 
1198  std::array<SF::Point, 3> tri_pts{base[i], base[j], base[N - 1]};
1199  std::array<double, 3> tri_vals{values[i], values[j], values[N - 1]};
1200 
1201  LocalSolver<3> solver(D, tri_pts, tri_vals);
1202  // min = std::min(min, solver.solve());
1203  SF_real tmp = solver.solve();
1204  if (min > tmp && tmp > T_R[indX] && compute_H(indX, tmp) > 0.0) {
1205  min = tmp;
1206  }
1207  }
1208  }
1209  }
1210 
1211  // 3) edges that include the target
1212  for (int i = 0; i < N - 1; i++) {
1213  std::array<int, 2> edge_ids{nodeIDs[i], nodeIDs[N - 1]};
1214  if (is_not_valid_update<2>(edge_ids, time)) continue;
1215 
1216  std::array<SF::Point, 2> edge_pts{base[i], base[N - 1]};
1217  std::array<double, 2> edge_vals{values[i], values[N - 1]};
1218 
1219  LocalSolver<2> solver(D, edge_pts, edge_vals);
1220  SF_real tmp = solver.solve();
1221  if (min > tmp && tmp > T_R[indX] && compute_H(indX, tmp) > 0.0) {
1222  min = tmp;
1223  }
1224  }
1225  }
1226 
1227  return min;
1228 }
1229 
1230 template <int N>
1231 bool eikonal_solver::is_not_valid_update(const std::array<int, N>& nodeIDs, double time)
1232 {
1233  const int target = nodeIDs[N - 1];
1234  const bool failedStim = (stim_status[target] == 2);
1235 
1236  for (int i = 0; i < N - 1; i++) {
1237  const int nb = nodeIDs[i];
1238 
1239  if (T_A[nb] < time) return true;
1240  if (T_A[nb] < T_R[nb]) return true;
1241  if (failedStim && stim_status[nb] == 1) return true;
1242  }
1243 
1244  return false;
1245 }
1246 
1247 SF_real eikonal_solver::compute_H(mesh_int_t& indX, SF_real& tmpTA)
1248 {
1249  // Apply a refractory delay if stimulus previously failed
1250  const double delay = (stim_status[indX] == 2) ? 5.0 : 0.0;
1251 
1252  // Early exit conditions
1253  if (tmpTA == inf ||
1254  (T_R[indX] + delay) == -400.0 ||
1255  std::fabs(tmpTA - TA_old[indX]) < param_globals::dream.fim.tol ||
1256  T_R[indX] > tmpTA) {
1257  return 1.0;
1258  }
1259 
1260  // Compute diastolic interval
1261  const double DI = tmpTA - (T_R[indX] + delay);
1262  D_I[indX] = DI;
1263 
1264  // CV restitution factor
1265  const double exponent = -(DI + kappa_cvrest[indX]) * denom_cvrest[indX];
1266  const double Factor_DI = 1.0 - rho_cvrest[indX] * std::exp(exponent);
1267 
1268  if (Factor_DI <= 0.0) {
1269  return 0.0;
1270  }
1271 
1272  // Apply restitution curve behavior
1273  return (DI <= theta_cvrest[indX]) ? -Factor_DI : Factor_DI;
1274 }
1275 
1276 SF_real eikonal_solver::compute_coherence(mesh_int_t& indX)
1277 {
1278  SF_real scaling = 1.0;
1279  SF_real p = T_A[indX]; // starting guess (previous arrival time)
1280  SF_real q = update(indX, scaling, true); // candidate new arrival time
1281 
1282  for (int niter = 0; niter < param_globals::dream.fim.max_coh; ++niter) {
1283  scaling = compute_H(indX, p); // restitution-based scaling for CV (can be negative mid-iteration)
1284  q = update(indX, std::fabs(scaling), true); // new candidate arrival time; use fabs() to allow intermediate negative scaling
1285 
1286  if (std::fabs(p - q) < param_globals::dream.fim.tol)
1287  break; // converged
1288 
1289  p = q; // continue iterating
1290  }
1291 
1292  // Only accept q if final state is physiologically valid
1293  return (compute_H(indX, q) < 0.0) ? T_A[indX] : q;
1294 }
1295 
1296 void eikonal_solver::compute_bc()
1297 {
1298  bool EmpList = sum(List) == 0;
1299 
1300  if (EmpList) {
1301  for (size_t j = Index_currStim; j < StimulusTimes.size(); j++) {
1302  if ((StimulusPoints[j] == -1) || (StimulusTimes[j] > StimulusTimes[Index_currStim])) {
1303  Index_currStim = j;
1304  break;
1305  }
1306 
1307  mesh_int_t indNode = StimulusPoints[j];
1308  SF_real TimeSt = StimulusTimes[j];
1309 
1310  bool cond2add = add_node_neighbor_to_list(T_R[indNode], T_A[indNode], TimeSt);
1311  if (cond2add && compute_H(indNode, TimeSt) > 0) {
1312  T_A[indNode] = TimeSt;
1313  stim_status[indNode] = 1;
1314  stats.bc_status = true;
1315 
1316  if (param_globals::dream.output.debugNode == indNode) {
1317  nodeData.T_A = TimeSt;
1319  }
1320 
1321  for (int ii = n2n_dsp[indNode]; ii < n2n_dsp[indNode + 1]; ii++) {
1322  mesh_int_t indXNB = n2n_connect[ii];
1323  if (List[indXNB] == 0) {
1324  List[indXNB] = 1;
1325  if (param_globals::dream.output.debugNode == indXNB) {
1327  nodeData.idXNB = indNode;
1328  nodeData.nbn_T_A = TimeSt;
1329  }
1330  }
1331  }
1332  }
1333 
1334  if (cond2add && compute_H(indNode, TimeSt) <= 0) {
1335  stim_status[indNode] = 2;
1336  }
1337  }
1338 
1339  // After NB of initial points are assigned remove initial points from the list if they were added in the previous loop.
1340  for (size_t j = 0; j < StimulusPoints.size(); j++) {
1341  if ((StimulusPoints[j] == -1) && (StimulusTimes[j] > StimulusTimes[Index_currStim])) continue;
1342  if (List[StimulusPoints[j]] == 1) {
1343  List[StimulusPoints[j]] = 0;
1344  if (param_globals::dream.output.debugNode == StimulusPoints[j]) nodeData.update_status(node_stats::out, node_stats::stim);
1345  }
1346  }
1347 
1349  } else {
1350  for (size_t i = Index_currStim; i < StimulusTimes.size(); i++) {
1351  mesh_int_t indNode = StimulusPoints[i];
1352  SF_real TimeSt = StimulusTimes[i];
1353 
1354  if (indNode == -1) continue;
1355 
1356  if (TimeSt <= actMAX) {
1357  bool cond2add = add_node_neighbor_to_list(T_R[indNode], T_A[indNode], TimeSt);
1358  if (cond2add && compute_H(indNode, TimeSt) > 0) {
1359  T_A[indNode] = TimeSt;
1360  stim_status[indNode] = 1;
1361  stats.bc_status = true;
1362 
1363  if (param_globals::dream.output.debugNode == indNode) {
1364  nodeData.T_A = TimeSt;
1366  }
1367 
1368  for (int ii = n2n_dsp[indNode]; ii < n2n_dsp[indNode + 1]; ii++) {
1369  mesh_int_t indXNB = n2n_connect[ii];
1370  if (List[indXNB] == 0) {
1371  List[indXNB] = 1;
1372  if (param_globals::dream.output.debugNode == indXNB) {
1374  nodeData.idXNB = indNode;
1375  nodeData.nbn_T_A = TimeSt;
1376  }
1377  }
1378  }
1379 
1380  for (size_t j = 0; j < StimulusPoints.size(); j++) {
1381  if ((StimulusPoints[j] == -1) && (StimulusTimes[j] > StimulusTimes[Index_currStim])) continue;
1382  if (List[StimulusPoints[j]] == 1) {
1383  List[StimulusPoints[j]] = 0;
1384  if (param_globals::dream.output.debugNode == StimulusPoints[j]) nodeData.update_status(node_stats::out, node_stats::stim);
1385  }
1386  }
1387  }
1388  if (cond2add && compute_H(indNode, TimeSt) <= 0) {
1389  stim_status[indNode] = 2;
1390  }
1391 
1392  } else {
1393  Index_currStim = i;
1394 
1395  break;
1396  }
1397  }
1398  }
1399 
1400 } // end compute stimulus
1401 
1403 {
1404  bool act_is_in_safety_window, empty_list_and_illegal_stimulus, empty_stimulus, stimulus_is_in_safety_window;
1405 
1406  act_is_in_safety_window = (actMIN > param_globals::dream.tau_s) && (actMIN - time > param_globals::dream.tau_s);
1407  empty_stimulus = (StimulusPoints[Index_currStim] == -1);
1408  stimulus_is_in_safety_window = (StimulusTimes[Index_currStim] - time > param_globals::dream.tau_s);
1409  empty_list_and_illegal_stimulus = (sum(List) == 0) && (empty_stimulus || stimulus_is_in_safety_window);
1410 
1411  return act_is_in_safety_window || empty_list_and_illegal_stimulus;
1412 }
1413 
1414 void eikonal_solver::update_Ta_in_active_list()
1415 {
1416  bool First_in_List = 1;
1417 
1418  for (size_t j = 0; j < List.size(); j++) {
1419  if (T_A[j] < 0) continue;
1420  if (List[j] == 1 && First_in_List) {
1421  actMIN = T_A[j];
1422  actMAX = T_A[j];
1423  First_in_List = 0;
1424  continue;
1425  }
1426  if (List[j] == 1) {
1427  if (T_A[j] < actMIN) actMIN = T_A[j];
1428  if (T_A[j] > actMAX) actMAX = T_A[j];
1429  }
1430  }
1431 }
1432 
1434 {
1435  for (size_t j = 0; j < List.size(); j++) {
1436  if (T_A[j] < user_globals::tm_manager->time) {
1437  TA_old[j] = T_A[j];
1438  stim_status[j] = 0;
1439  nReadded2List[j] = 0;
1440 
1441  if (List[j] == 1) List[j] = 0;
1442  }
1443  }
1444 }
1445 
1446 bool eikonal_solver::add_node_neighbor_to_list(SF_real& RT, SF_real& oldTA, SF_real& newTA)
1447 {
1448  // Condition (A): RT < newTA < oldTA
1449  // Condition (B): oldTA < RT < newTA
1450  return ((RT < newTA) && (newTA < oldTA)) ||
1451  ((oldTA < RT) && (RT < newTA));
1452 }
1453 
1455 {
1456  SF_real* old_ptr = Vmv_old.ptr();
1457  SF_real* new_ptr = Vmv.ptr();
1458 
1459  const double thresh = param_globals::dream.repol_time_thresh;
1460 
1461  for (int ind_nodes = 0; ind_nodes < num_pts; ind_nodes++) {
1462  if (old_ptr[ind_nodes] >= thresh && new_ptr[ind_nodes] < thresh) {
1463  T_R[ind_nodes] = time;
1464  }
1465  }
1466 
1467  Vmv_old.release_ptr(old_ptr);
1468  Vmv.release_ptr(new_ptr);
1469 }
1470 
1472 {
1473  double t1, t0;
1474  get_time(t0);
1475 
1476  #ifdef _OPENMP
1477  int max_threads = omp_get_max_threads();
1478  omp_set_num_threads(1); // with the current setup of this function, it is better to run it in serial.
1479  #endif
1480 
1481  limpet::MULTI_IF* miif = ion.miif;
1482 
1483  const double time = user_globals::tm_manager->time,
1485 
1486  for (int i = 0; i < miif->N_IIF; i++) {
1487  if (!miif->N_Nodes[i]) continue;
1488  limpet::IonIfBase* pIF = ion.miif->IIF[i];
1489  limpet::IonIfBase* IIF_old = pIF->get_type().make_ion_if(pIF->get_target(),
1490  pIF->get_num_node(),
1491  ion.miif->plugtypes[i]);
1492  IIF_old->copy_SVs_from(*pIF, false);
1493 
1494  int current = 0;
1495  do {
1496  int ind_gb = (miif->NodeLists[i][current]);
1497  // Global index of current node;
1498  // save states at the moment
1499  double Vm_old = miif->ldata[i][limpet::Vm][current];
1500  double Iion_old = miif->ldata[i][limpet::Iion][current];
1501  double prev_R = T_R[ind_gb];
1502 
1503  bool cond_2upd = T_R[ind_gb] <= T_A[ind_gb] && T_A[ind_gb] != inf && Vm_old > param_globals::dream.repol_time_thresh;
1504 
1505  if (cond_2upd) {
1506  double elapsed_time = 0;
1507  double prev_Vm = miif->ldata[i][limpet::Vm][current];
1508 
1509  do {
1510  elapsed_time += dt;
1511  pIF->compute(current, current + 1, miif->ldata[i]);
1512  miif->ldata[i][limpet::Vm][current] -= miif->ldata[i][limpet::Iion][current] * param_globals::dt;
1513  if (miif->ldata[i][limpet::Vm][current] < param_globals::dream.repol_time_thresh) {
1514  T_R[ind_gb] = time + elapsed_time;
1515  break;
1516  }
1517  } while (elapsed_time < 600);
1518 
1519  // Restore old states
1520  miif->ldata[i][limpet::Vm][current] = Vm_old;
1521  miif->ldata[i][limpet::Iion][current] = Iion_old;
1522  }
1523  current++;
1524  } while (current < miif->N_Nodes[i]);
1525  pIF->copy_SVs_from(*IIF_old, false);
1526  }
1527 
1528  #ifdef _OPENMP
1529  omp_set_num_threads(max_threads); // restore max threads for parabolic solver
1530  #endif
1531 
1532  // treat solver statistics
1533  auto dur = timing(t1, t0);
1534  stats.slvtime_D += dur;
1535 }
1536 
1538 {
1539  double t0, t1;
1540  get_time(t0);
1541 
1542  SF_real* c = Idiff->ptr();
1543  SF_real* v = vm.ptr();
1544 
1545  const sf_mesh& mesh = get_mesh(intra_elec_msh);
1546  const SF::vector<mesh_int_t>& alg_nod = mesh.pl.algebraic_nodes();
1547  int rank = get_rank();
1548 
1549  for (size_t j = 0; j < alg_nod.size(); j++) {
1550  mesh_int_t loc_nodal_idx = alg_nod[j];
1551  mesh_int_t loc_petsc_idx = local_nodal_to_local_petsc(mesh, rank, loc_nodal_idx);
1552 
1553  double TA = T_A[loc_nodal_idx];
1554  if (!(time >= TA && (TA + 5) >= time && List[loc_nodal_idx] == 0))
1555  continue;
1556 
1557  double dT = time - TA;
1558  auto& node = diff_cur[loc_nodal_idx];
1559 
1560  switch (node.model) {
1561  case GAUSS: {
1562  double term1 = (dT - node.beta_1) / node.gamma_1;
1563  double term2 = (dT - node.beta_2) / node.gamma_2;
1564  double term3 = (dT - node.beta_3) / node.gamma_3;
1565  c[loc_petsc_idx] = node.alpha_1 * exp(-term1 * term1) + node.alpha_2 * exp(-term2 * term2) + node.alpha_3 * exp(-term3 * term3);
1566  break;
1567  }
1568  default: {
1569  double e_on = (dT >= 0.0) ? 1.0 : 0.0;
1570  double e_off = (v[loc_petsc_idx] < node.V_th) ? 1.0 : 0.0;
1571  c[loc_petsc_idx] = node.A_F / node.tau_F * exp(dT / node.tau_F) * e_on * e_off;
1572  break;
1573  }
1574  }
1575  }
1576 
1577  Idiff->release_ptr(c);
1578  vm.release_ptr(v);
1579 
1580  stats.slvtime_B += timing(t1, t0);
1581 }
1582 
1583 void eikonal_solver::save_eikonal_state(const char* tsav_ext)
1584 {
1585  if (get_rank() == 0) {
1586  FILE* file_writestate;
1587  char buffer_writestate[1024];
1588  snprintf(buffer_writestate, sizeof buffer_writestate, "%s.%s.roe.dat", param_globals::write_statef, tsav_ext);
1589  file_writestate = fopen(buffer_writestate, "w");
1590  for (size_t jjj = 0; jjj < List.size(); jjj++) {
1591  fprintf(file_writestate, "%d %d %d %f %f %f %f \n", List[jjj], num_changes[jjj], nReadded2List[jjj], T_A[jjj], T_R[jjj], TA_old[jjj], D_I[jjj]);
1592  }
1593 
1594  fclose(file_writestate);
1595  }
1596 }
1597 
1598 void eikonal_solver::load_state_file()
1599 {
1600  set_dir(INPUT);
1601  FILE* file_startstate;
1602  char buffer_startstate[strlen(param_globals::start_statef) + 10];
1603 
1604  snprintf(buffer_startstate, sizeof buffer_startstate, "%s.dat", param_globals::start_statef);
1605  file_startstate = fopen(buffer_startstate, "r");
1606 
1607  if (file_startstate == NULL) {
1608  log_msg(NULL, 5, 0, "Not able to open state file: %s", buffer_startstate);
1609  } else if (param_globals::output_level) {
1610  log_msg(NULL, 0, 0, "Open state file for eikonal model: %s", buffer_startstate);
1611  }
1612 
1613  int ListVal, numChangesVal, numChanges2Val;
1614  float T_AVal, T_RVal, TA_oldVal, D_IVal, PCLVal;
1615  int index = 0;
1616 
1617  while (fscanf(file_startstate, "%d %d %d %f %f %f %f", &ListVal, &numChangesVal, &numChanges2Val, &T_AVal, &T_RVal, &TA_oldVal, &D_IVal) == 7) {
1618  List[index] = ListVal;
1619  num_changes[index] = numChangesVal;
1620  nReadded2List[index] = numChanges2Val;
1621  T_A[index] = T_AVal;
1622  T_R[index] = T_RVal;
1623  TA_old[index] = TA_oldVal;
1624  D_I[index] = D_IVal;
1625  ++index;
1626  }
1627  if (param_globals::output_level) log_msg(NULL, 0, 0, "Number of nodes in active list: %i", sum(List));
1628 
1629  fclose(file_startstate);
1630 }
1631 
1633 {
1634  logger = f_open(filename, "w");
1635 
1636  const char* h1 = " ------ ---------- ---------- ------- ------- | List logic ----- -------- | Neighbor node ---- |";
1637  const char* h2 = " cycle AT old AT RT DI | Status Entry Exit | ID AT |";
1638 
1639  if (logger == NULL)
1640  log_msg(NULL, 3, 0, "%s error: Could not open file %s in %s. Turning off logging.\n",
1641  __func__, filename);
1642  else {
1643  log_msg(logger, 0, 0, "%s", h1);
1644  log_msg(logger, 0, 0, "%s", h2);
1645  }
1646 }
1647 
1648 void node_stats::log_stats(double time, bool cflg)
1649 {
1650  if (!this->logger) return;
1651 
1652  char abuf[256];
1653  char bbuf[256];
1654  char cbuf[256];
1655 
1656  // create nicer output in logger for inf, -inf, nan
1657  std::ostringstream oss_TA, oss_TA_, oss_TR, oss_DI, oss_nbnTA;
1658  oss_TA << this->T_A;
1659  oss_TA_ << this->T_A_;
1660  oss_TR << this->T_R;
1661  oss_DI << this->D_I;
1662  oss_nbnTA << this->nbn_T_A;
1663 
1664  if (this->idXNB == std::numeric_limits<Int>::min()) {
1665  snprintf(cbuf, sizeof cbuf, "%7s %10s", "-", "-");
1666  } else {
1667  snprintf(cbuf, sizeof cbuf, "%7d %10s", this->idXNB, oss_nbnTA.str().c_str());
1668  }
1669 
1670  snprintf(abuf, sizeof abuf, "%6d %10s %10s %7s %7s", this->cycle, oss_TA.str().c_str(), oss_TA_.str().c_str(), oss_TR.str().c_str(), oss_DI.str().c_str());
1671  snprintf(bbuf, sizeof bbuf, "%7s %8s %8s", this->status, this->reasonIn, this->reasonOut);
1672 
1673  unsigned char flag = cflg ? ECHO : 0;
1674  log_msg(this->logger, 0, flag | FLUSH | NONL, "%9.3f %s | %s | %s |\n", time, abuf, bbuf, cbuf);
1675 
1676  this->reasonIn = "-";
1677  this->reasonOut = "-";
1678  this->T_A_ = this->T_A;
1679  this->T_A = std::numeric_limits<double>::quiet_NaN();
1680  this->T_R = std::numeric_limits<double>::quiet_NaN();
1681  this->D_I = std::numeric_limits<double>::quiet_NaN();
1683  this->nbn_T_A = std::numeric_limits<double>::quiet_NaN();
1684  this->cycle++;
1685 }
1686 
1688 {
1689  // directly convert enums to strings for logging
1690  const char* stat_str;
1691  const char* reas_str;
1692  switch (s) {
1693  case node_stats::in: stat_str = "in"; break;
1694  case node_stats::out: stat_str = "out"; break;
1695  }
1696 
1697  switch (r) {
1698  case node_stats::none: reas_str = "-"; break;
1699  case node_stats::nbn: reas_str = "nbn"; break;
1700  case node_stats::conv: reas_str = "conv"; break;
1701  case node_stats::stim: reas_str = "stim"; break;
1702  }
1703 
1704  this->status = stat_str;
1705  if (s == in) {
1706  this->reasonIn = reas_str;
1707  } else {
1708  this->reasonOut = reas_str;
1709  }
1710 }
1711 
1712 template<int MESH_SIZE>
1714  const SF::Point& x1 = points[0];
1715  const SF::Point& x2 = points[1];
1716  const SF::Point& x3 = points[2];
1717  const SF::Point& x4 = points[3];
1718 
1719  const double& u1 = values[0];
1720  const double& u2 = values[1];
1721  const double& u3 = values[2];
1722 
1723  if constexpr (MESH_SIZE == 2) {
1724  return tsitsiklis_update_line({x1,x2}, D, u1);
1725 
1726  } else if constexpr (MESH_SIZE == 3) {
1727  return tsitsiklis_update_triangle({x1,x2,x3}, D, {u1,u2});
1728 
1729  } else if constexpr (MESH_SIZE == 4) {
1730  double u_tet = tsitsiklis_update_tetra({x1,x2,x3,x4}, D, {u1,u2,u3});
1731  if (isnan(u_tet)) {
1732  u_tet = std::numeric_limits<double>::infinity();
1733  }
1734  // face calculations (contains edge update as fallback)
1735  double u_face1 = tsitsiklis_update_triangle({x1, x2, x4}, D, {u1, u2});
1736  double u_face2 = tsitsiklis_update_triangle({x1, x3, x4}, D, {u1, u3});
1737  double u_face3 = tsitsiklis_update_triangle({x2, x3, x4}, D, {u2, u3});
1738 
1739  double u_tri = std::min({u_face1, u_face2, u_face3});
1740  return std::min(u_tet, u_tri);
1741  }
1742 }
1743 
1744 template <int MESH_SIZE>
1745 double LocalSolver<MESH_SIZE>::tsitsiklis_update_line(const std::array<SF::Point, 2>& base,
1746  const SF::dmat<double>& D,
1747  const double& value)
1748 {
1749  // Compute the difference vector: a1 = x2 - x1
1750  SF::Point a1 = base[1] - base[0];
1751  // Return updated value at x2: u1 + norm
1752  return value + std::sqrt(SF::inner_prod(a1, D * a1));
1753 }
1754 
1755 template <int MESH_SIZE>
1756 double LocalSolver<MESH_SIZE>::tsitsiklis_update_triangle(const std::array<SF::Point, 3>& base,
1757  const SF::dmat<double>& D,
1758  const std::array<double, 2>& values)
1759 {
1760  const SF::Point& x1 = base[0];
1761  const SF::Point& x2 = base[1];
1762  const SF::Point& x3 = base[2];
1763  const double& u1 = values[0];
1764  const double& u2 = values[1];
1765  double result = std::numeric_limits<double>::infinity();
1766 
1767  SF::Point z1 = x1 - x2;
1768  SF::Point z2 = x2 - x3;
1769  double k = u1 - u2;
1770 
1771  // Squared Mahalanobis norms
1772  const SF::Point Dz1 = D * z1;
1773  const SF::Point Dz2 = D * z2;
1774  double p11 = SF::inner_prod(z1, Dz1);
1775  double p12 = SF::inner_prod(z1, Dz2);
1776  double p22 = SF::inner_prod(z2, Dz2);
1777 
1778  double denominator = p11 - k * k;
1779  double sqrt_val = (p11 * p22 - p12 * p12) / denominator;
1780 
1781  if (denominator > 1e-12) { // avoid dividing by zero or near zero -> k very small, likely due to collapsed triangle/coinciding points
1782  const double sqrt_val = (p11 * p22 - p12 * p12) / denominator;
1783  if (sqrt_val >= 0.0) { // only real solutions are considered
1784  const double rhs = k * std::sqrt(sqrt_val);
1785  double alpha1 = -(p12 + rhs) / p11;
1786  double alpha2 = -(p12 - rhs) / p11;
1787 
1788  alpha1 = std::clamp(alpha1, 0.0, 1.0);
1789  alpha2 = std::clamp(alpha2, 0.0, 1.0);
1790 
1791  for (double alpha : {alpha1, alpha2}) {
1792  SF::Point x_interp = x1 * alpha + x2 * (1.0 - alpha);
1793  SF::Point dist = x3 - x_interp;
1794  double norm_D = std::sqrt(SF::inner_prod(dist, D * dist));
1795  double u3 = alpha * u1 + (1.0 - alpha) * u2 + norm_D;
1796  result = std::min(result, u3);
1797  }
1798  }
1799  }
1800 
1801  // Fallback: point-based update if square root was invalid or denominator is too small
1802  double u_edge1 = tsitsiklis_update_line({x1,x3}, D, u1);
1803  double u_edge2 = tsitsiklis_update_line({x2,x3}, D, u2);
1804 
1805  return std::min({result, u_edge1, u_edge2});
1806 }
1807 
1808 template <int MESH_SIZE>
1809 double LocalSolver<MESH_SIZE>::tsitsiklis_update_tetra(const std::array<SF::Point, 4>& base,
1810  const SF::dmat<double>& D,
1811  const std::array<double, 3>& values)
1812 {
1813  const SF::Point& x1 = base[0];
1814  const SF::Point& x2 = base[1];
1815  const SF::Point& x3 = base[2];
1816  const SF::Point& x4 = base[3];
1817 
1818  const double& u1 = values[0];
1819  const double& u2 = values[1];
1820  const double& u3 = values[2];
1821 
1822  // edge vectors
1823  const SF::Point y1 = x3 - x1;
1824  const SF::Point y2 = x3 - x2;
1825  const SF::Point y3 = x4 - x3;
1826 
1827  const double k1 = u1 - u3;
1828  const double k2 = u2 - u3;
1829 
1830  // Matrix products (squared norms and dot products under metric D)
1831  const SF::Point Dy1 = D * y1;
1832  const SF::Point Dy2 = D * y2;
1833  const SF::Point Dy3 = D * y3;
1834  const double r11 = SF::inner_prod(y1, Dy1);
1835  const double r12 = SF::inner_prod(y1, Dy2);
1836  const double r13 = SF::inner_prod(y1, Dy3);
1837  const double r21 = r12;
1838  const double r22 = SF::inner_prod(y2, Dy2);
1839  const double r23 = SF::inner_prod(y2, Dy3);
1840  const double r31 = r13;
1841  const double r32 = r23;
1842 
1843  const double A1 = k2 * r11 - k1 * r12;
1844  const double A2 = k2 * r21 - k1 * r22;
1845  const double B = k2 * r31 - k1 * r32;
1846  const double k = k1 - (A1 / A2) * k2;
1847  const SF::Point z1 = y1 - (A1 / A2) * y2;
1848  const SF::Point z2 = y3 - (B / A2) * y2;
1849 
1850  // compute quadratic equation
1851  const SF::Point Dz1 = D * z1;
1852  const SF::Point Dz2 = D * z2;
1853  const double p11 = SF::inner_prod(z1, Dz1);
1854  const double p12 = SF::inner_prod(z1, Dz2);
1855  const double p22 = SF::inner_prod(z2, Dz2);
1856  const double denominator = p11 - k*k;
1857  const double sqrt_val = (p11 * p22 - (p12 * p12)) / denominator;
1858  const double rhs = k * std::sqrt(sqrt_val);
1859 
1860  double alpha1 = -(p12 + rhs) / p11;
1861  double alpha2 = -(B + alpha1 * A1) / A2;
1862 
1863  // handle degenerate cases
1864  const double EPS = 1e-16;
1865  if ((std::abs(A1) < EPS) && (std::abs(A2) < EPS)) {
1866  alpha1 = (r12 * r23 - r13 * r22) / (r11 * r22 - (r12 * r12));
1867  alpha2 = (r12 * r13 - r11 * r23) / (r11 * r22 - (r12 * r12));
1868  } else if ((std::abs(A1) < EPS) && (std::abs(A2) > EPS)) {
1869  alpha1 = 0;
1870  alpha2 = -B / A2;
1871  } else if ((std::abs(A1) > EPS) && (std::abs(A2) < EPS)) {
1872  alpha1 = -B / A1;
1873  alpha2 = 0;
1874  }
1875 
1876  double alpha3 = 1 - alpha1 - alpha2;
1877  const SF::Point dist = x4 - (alpha1 * x1 + alpha2 * x2 + alpha3 * x3);
1878  // barycentric coordinate should be inside the tetrahedron
1879  if (alpha1 < -EPS || alpha2 < -EPS || alpha3 < -EPS ||
1880  alpha1 > 1.0+EPS || alpha2 > 1.0+EPS || alpha3 > 1.0+EPS) {
1881  return std::numeric_limits<double>::infinity();
1882  } else {
1883  return alpha1 * u1 + alpha2 * u2 + alpha3 * u3 + std::sqrt(SF::inner_prod(dist, D * dist));
1884  }
1885 }
1886 
1887 } // namespace opencarp
void output(vector< int > nodes, IGBheader *h, char *ofname, enum_format of, int t0, int t1, int stride, bool explode, float scale)
Definition: IGBextract.cc:208
int mesh_int_t
Definition: SF_container.h:46
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
Basic utility structs and functions, mostly IO related.
#define FLUSH
Definition: basics.h:311
#define ECHO
Definition: basics.h:308
#define NONL
Definition: basics.h:312
virtual void write(const char *filename) const =0
virtual S * ptr()=0
virtual void release_ptr(S *&p)=0
virtual void add_scaled(const abstract_vector< T, S > &vec, S k)=0
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:429
vector< T > dsp
connectivity starting index of each element
Definition: SF_container.h:416
vector< S > she
sheet direction
Definition: SF_container.h:421
vector< elem_t > type
element type
Definition: SF_container.h:418
vector< T > con
Definition: SF_container.h:412
size_t l_numpts
local number of points
Definition: SF_container.h:401
Container for a PETSc VecScatter.
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
const T * end() const
Pointer to the vector's end.
Definition: SF_vector.h:128
void assign(InputIterator s, InputIterator e)
Assign a memory range.
Definition: SF_vector.h:161
void reserve(size_t n)
Definition: SF_vector.h:241
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
Represents the ionic model and plug-in (IMP) data structure.
Definition: ION_IF.h:139
const IonType & get_type() const
Gets this IMP's model type.
Definition: ION_IF.cc:144
virtual void copy_SVs_from(IonIfBase &other, bool alloc)=0
Copies the state variables of an IMP.
Target get_target() const
Definition: ION_IF.h:339
int get_num_node() const
Gets the number of nodes handled by this IMP.
Definition: ION_IF.cc:148
void compute(int start, int end, GlobalData_t **data)
Perform ionic model computation for 1 time step.
Definition: ION_IF.cc:270
virtual IonIfBase * make_ion_if(Target target, int num_node, const std::vector< std::reference_wrapper< IonType >> &plugins) const =0
Generate an IonIf object from this type.
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
std::vector< IonTypeList > plugtypes
plugins types for each region
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: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
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
const char * name
The name of the physic, each physic should have one.
Definition: physics_types.h:62
std::string timer_unit(const int timer_id)
figure out units of a signal linked to a given timer
SF::vector< stimulus > stimuli
the electrical stimuli
parabolic_solver parab_solver
Solver for the parabolic bidomain equation.
MaterialType mtype[2]
the material types of intra_grid and extra_grid grids.
void destroy()
Currently we only need to close the file logger.
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
sf_vec * phie_dummy
no elliptic solver needed, but we need a dummy for phie to use parabolic solver
eikonal_solver eik_solver
Solver for the eikonal equation.
LAT_detector lat
the activation time detector
gvec_data gvec
datastruct holding global IMP state variable output
generic_timing_stats IO_stats
grid_t
An electrics grid identifier to distinguish between intra and extra grids.
igb_output_manager output_manager_cycle
void initialize()
Initialize the Eikonal class.
igb_output_manager output_manager_time
class handling the igb output
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
int check_quiescence(double tm, double dt)
check for quiescence
Definition: electrics.cc:1725
void output_initial_activations()
output one nodal vector of initial activation time
Definition: electrics.cc:1834
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:1550
int check_acts(double tm)
check activations at sim time tm
Definition: electrics.cc:1670
void FIM()
Standard fast iterative method to solve eikonal equation with active list approach.
void update_repolarization_times_from_rd(sf_vec &Vmv, sf_vec &Vmv_old, double time)
Updates node repolarization times based on transmembrane voltage crossing.
SF::vector< SF_real > T_R
void init_imp_region_properties()
Initializes diffusion current models and CV restitution parameters per mesh node.
SF::vector< mesh_int_t > n2e_dsp
SF::vector< SF_real > T_A
SF::vector< SF_real > TA_old
SF::vector< diffusion_current > diff_cur
void init()
Initialize vectors and variables in the eikonal_solver class.
SF::vector< mesh_int_t > e2n_con
bool determine_model_to_run(double &time)
Determine the next model to run in the alternation between RD and Eikonal.
SF::vector< mesh_int_t > stim_status
SF::vector< mesh_int_t > elem_start
SF::vector< mesh_int_t > StimulusPoints
SF::vector< SF_real > rho_cvrest
std::vector< double > CV_L
void save_eikonal_state(const char *tsav_ext)
Save the current state of variables related to the Eikonal simulation to a file to initialize a futur...
SF::vector< SF_real > denom_cvrest
void compute_diffusion_current(const double &time, sf_vec &vm)
Computes the stimulus-driven diffusion current at mesh nodes.
std::vector< mesh_int_t > n2n_connect
void set_stimuli(SF::vector< stimulus > &stimuli)
Simple setter for stimulus vector.
SF::vector< mesh_int_t > e2n_cnt
SF::vector< mesh_int_t > n2e_con
void clean_list()
Clean the list of nodes by resetting their status and tracking changes based on the time step of the ...
SF::vector< SF_real > D_I
SF::vector< mesh_int_t > nReadded2List
SF::vector< SF_real > theta_cvrest
SF::vector< mesh_int_t > num_changes
std::vector< mesh_int_t > n2n_dsp
SF::vector< SF_real > StimulusTimes
SF::vector< SF_real > kappa_cvrest
SF::vector< mesh_int_t > n2e_cnt
void cycFIM()
Implementation of the cyclical fast iterative method used in step A of the DREAM model.
void update_repolarization_times(const Ionics &ion)
Estimates initial repolarization times (T_R) in Step D of DREAM.
eikonal_solver_stats stats
std::vector< SF::dmat< double > > S
void write_data()
write registered data to disk
Definition: sim_utils.cc:1494
void close_files_and_cleanup()
close file descriptors
Definition: sim_utils.cc:1541
void register_output(sf_vec *inp_data, const mesh_t inp_meshid, const int dpn, const char *name, const char *units, const SF::vector< mesh_int_t > *idx=NULL, bool elem_data=false)
Register a data vector for output.
Definition: sim_utils.cc:1461
sf_mat * rhs_parab
rhs matrix to solve parabolic
Definition: electrics.h:119
lin_solver_stats stats
Definition: electrics.h:129
void rebuild_matrices(MaterialType *mtype, limpet::MULTI_IF &miif, FILE_SPEC logger)
Definition: electrics.cc:1202
void solve(sf_vec &phie_i)
Definition: electrics.cc:1315
sf_mat * mass_i
lumped for parabolic problem
Definition: electrics.h:118
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 * old_vm
older Vm needed for 2nd order dT
Definition: electrics.h:110
int write_trace()
write traces to file
Definition: signals.h:701
stim_t type
type of stimulus
Definition: stimulate.h:137
int npls
number of stimulus pulses
Definition: stimulate.h:120
double pcl
pacing cycle length
Definition: stimulate.h:121
double start
start time of protocol
Definition: stimulate.h:119
sig::time_trace wave
wave form of stimulus pulse
Definition: stimulate.h:97
stim_protocol ptcl
applied stimulation protocol used
Definition: stimulate.h:168
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
bool is_active() const
Return whether stim is active.
Definition: stimulate.cc:190
void setup(int idx)
Setup from a param stimulus index.
Definition: stimulate.cc:166
void dump_vtx_file(int idx)
Export the vertices to vtx file.
Definition: stimulate.cc:457
stim_physics phys
physics of stimulus
Definition: stimulate.h:169
std::string name
label stimulus
Definition: stimulate.h:165
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
double time
current time
Definition: timer_utils.h:76
Diffusion Reaction Eikonal Alternant Model (DREAM) based on the electrics physics class.
void transpose_connectivity(const vector< T > &a_cnt, const vector< T > &a_con, vector< T > &b_cnt, vector< T > &b_con)
Transpose CRS matrix graph A into B.
double inner_prod(const Point &a, const Point &b)
Definition: SF_container.h:90
T sum(const vector< T > &vec)
Compute sum of a vector's entries.
Definition: SF_vector.h:340
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
void outer_prod(const Point &a, const Point &b, const double s, double *buff, const bool add=false)
Definition: SF_container.h:95
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 init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:99
V clamp(const V val, const W start, const W end)
Clamp a value into an interval [start, end].
Definition: kdpart.hpp:132
void cnt_from_dsp(const std::vector< T > &dsp, std::vector< T > &cnt)
Compute counts from displacements.
Definition: kdpart.hpp:149
void dsp_from_cnt(const std::vector< T > &cnt, std::vector< T > &dsp)
Compute displacements from counts.
Definition: kdpart.hpp:140
constexpr T min(T a, T b)
Definition: ion_type.h:33
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
int stimidx_from_timeridx(const SF::vector< stimulus > &stimuli, const int timer_id)
determine link between timer and stimulus
Definition: electrics.cc:800
@ iotm_chkpt_list
Definition: timer_utils.h:44
@ iotm_console
Definition: timer_utils.h:44
@ iotm_trace
Definition: timer_utils.h:44
@ iotm_chkpt_intv
Definition: timer_utils.h:44
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.
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
SF::scattering * get_permutation(const int mesh_id, const int perm_id, const int dpn)
Get the PETSC to canonical permutation scattering for a given mesh and number of dpn.
SF::meshdata< mesh_int_t, mesh_real_t > sf_mesh
Definition: sf_interface.h:48
void apply_stim_to_vector(const stimulus &s, sf_vec &vec, bool add)
Definition: electrics.cc:471
int set_dir(IO_t dest)
Definition: sim_utils.cc:490
vec3< V > cross(const vec3< V > &a, const vec3< V > &b)
Definition: vect.h:144
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:276
V dist(const vec3< V > &p1, const vec3< V > &p2)
Definition: vect.h:114
@ Vm_clmp
Definition: stimulate.h:78
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:363
void init_stim_info(void)
uses potential for stimulation
Definition: stimulate.cc:49
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:135
SF::scattering * register_permutation(const int mesh_id, const int perm_id, const int dpn)
Register a permutation between two orderings for a mesh.
void setup_dataout(const int dataout, std::string dataout_vtx, mesh_t grid, SF::vector< mesh_int_t > *&restr, bool async)
Definition: electrics.cc:631
@ OUTPUT
Definition: sim_utils.h:53
void init_sv_gvec(gvec_data &GVs, limpet::MULTI_IF *miif, sf_vec &tmpl, igb_output_manager &output_manager)
Definition: ionics.cc:567
void assemble_sv_gvec(gvec_data &gvecs, limpet::MULTI_IF *miif)
Definition: ionics.cc:638
char * dupstr(const char *old_str)
Definition: basics.cc:44
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
@ intra_elec_msh
Definition: sf_interface.h:60
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:63
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:290
const char * get_tsav_ext(double time)
Definition: electrics.cc:884
V timing(V &t2, const V &t1)
Definition: basics.h:448
void f_close(FILE_SPEC &f)
Close a FILE_SPEC.
Definition: basics.cc:162
@ ElecMat
Definition: fem_types.h:39
#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 ELEM_PETSC_TO_CANONICAL
Permute algebraic element data from PETSC to canonical ordering.
Definition: sf_interface.h:81
#define EXP_POSTPROCESS
Definition: sim_utils.h:162
Electrical stimulation functions.
Point and vector struct.
Definition: SF_container.h:65
double y
Definition: SF_container.h:67
double z
Definition: SF_container.h:68
double x
Definition: SF_container.h:66
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
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
double slvtime_A
total time in Step A
Definition: timers.h:45
void log_stats(double time, bool cflg)
Definition: timers.cc:130
double minAT
minimum activation time in current solve
Definition: timers.h:40
double maxAT
maximum activation time in current solve
Definition: timers.h:42
void init_logger(const char *filename)
Definition: timers.cc:114
int activeList
number of nodes currently in list
Definition: timers.h:38
void update_iter(const int curiter)
Definition: timers.cc:162
double slvtime_B
total time in Step B
Definition: timers.h:47
double slvtime_D
total time in Step D
Definition: timers.h:49
bool bc_status
boundary conditions were applied?
Definition: timers.h:53
void update_cli(double time, bool cflg)
Definition: timers.cc:168
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 log_stats(double tm, bool cflg)
Definition: timers.cc:27
const char * reasonOut
reason for list entry
SF_real T_R
repolarization time
void init_logger(const char *filename)
void log_stats(double tm, bool cflg)
SF_real D_I
diastolic interval
mesh_int_t idXNB
neighboring node index responsible for list entry
const char * reasonIn
reason for list entry
SF_real T_A_
previous activation time
SF_real T_A
current activation time
SF_real nbn_T_A
activation time of neighboring node
mesh_int_t cycle
DREAM cycle.
mesh_int_t idX
node index
void update_status(enum status s, enum reason r)