openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
sim_utils.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 "basics.h"
28 #include "sim_utils.h"
29 #include "fem.h"
30 #include "physics.h"
31 #include "async_io.h"
32 #include "SF_init.h"
33 
34 #include <libgen.h>
35 #include <fstream>
36 #include <iomanip>
37 
38 
39 namespace opencarp {
40 
41 static char input_dir[1024], // directory from which to read input
42  output_dir[1024], // directory to which to write results
43  postproc_dir[1024], // postprocessing directory
44  current_dir[1024]; // current directory
45 
46 void parse_params_cpy(int argc, char** argv)
47 {
48  // copy the command line parameters that are for carp (i.e. before the first "+")
49  int cpy_argc = argc;
50  SF::vector<char*> cpy_argv(cpy_argc, NULL);
51 
52  for(int i=0; i < cpy_argc; i++) {
53  if(strcmp(argv[i], "+") != 0) {
54  cpy_argv[i] = dupstr(argv[i]);
55  }
56  else {
57  cpy_argc = i;
58  break;
59  }
60  }
61 
62  // launch param on the carp command line parameters
63  int param_argc = cpy_argc;
64  int status;
65 
66  do {
67  status = param(PARAMETERS, &param_argc, cpy_argv.data());
68  if ( status==PrMERROR||status==PrMFATAL )
69  fprintf( stderr, "\n*** Error reading parameters\n\n");
70  else if ( status == PrMQUIT )
71  fprintf( stderr, "\n*** Quitting by user's request\n\n");
72  } while (status == PrMERROR);
73 
74  // if --output-setup then loop over PARAMETERS
75  // output value of each PARAMETER
76  if (param_globals::output_setup)
77  param_output(PARAMETERS, 1);
78 
79  // free the parameters
80  for(int i=0; i<param_argc; i++)
81  free(cpy_argv[i]);
82 
83  // exit on error
84  if (status & PrMERROR) exit(status);
85 
86 }
87 
88 
90 {
91  // here all the physics can be registered to the physics registry
92  // then they should be processed automatically
95 
98 
101 }
102 
104 {
105  log_msg(0,0,0, "\n *** Initializing physics ***\n");
106 
107  //load in the external imp modules
108 #ifdef HAVE_DLOPEN
109  for (int ii = 0; ii < param_globals::num_external_imp; ii++) {
110  int loading_succeeded = limpet::load_ionic_module(param_globals::external_imp[ii]);
111  assert(loading_succeeded);
112  }
113 #else
114  if(param_globals::num_external_imp)
115  log_msg(NULL, 4, ECHO,"Loading of external LIMPET modules not enabled.\n"
116  "Recompile with DLOPEN set.\n" );
117 #endif
118 
119  // init physics via Basic_physic interface
120  for(auto it : user_globals::physics_reg) {
121  Basic_physic* p = it.second;
122  log_msg(NULL, 0, 0, "Initializing %s ..", p->name);
123  p->initialize();
124  }
125 }
126 
128 {
129  log_msg(0,0,0, "\n *** Destroying physics ***\n");
130 
131  for(auto it : user_globals::physics_reg) {
132  Basic_physic* p = it.second;
133  log_msg(NULL, 0, 0, "Destroying %s ..", p->name);
134  p->destroy();
135  }
136 }
137 
138 // ignore_extracellular stim must be moved to stimulate.cc to be able
139 // to use all defined set operations instead of defines
140 
154 void ignore_extracellular_stim(Stimulus *st, int ns, int ignore)
155 {
156  // needs to be switch to stim enum types defined in stimulate.h
157  for ( int i=0; i<ns; i++ ) {
158  int turn_off = 0;
159  turn_off += (st[i].stimtype == Extracellular_Ground) && (ignore & NO_EXTRA_GND);
160  turn_off += (IsExtraV(st[i])) && (ignore & NO_EXTRA_V);
161  turn_off += (st[i].stimtype==Extracellular_I) && (ignore & NO_EXTRA_I);
162 
163  if (turn_off) {
164  st[i].stimtype = Ignore_Stim;
165  log_msg( NULL, 1, 0, "Extracellular stimulus %d ignored for monodomain", i );
166  } else if ( st[i].stimtype==Intracellular_I ) {
167  st[i].stimtype = Transmembrane_I;
168  log_msg( NULL, 1, 0, "Intracellular stimulus %d converted to transmembrane", i );
169  }
170  }
171 }
172 
180 int set_ignore_flags( int mode )
181 {
182  if(mode==MONODOMAIN)
183  return STM_IGNORE_MONODOMAIN;
184  if(mode==BIDOMAIN)
185  return STM_IGNORE_BIDOMAIN;
186  if(mode==PSEUDO_BIDM)
187  return STM_IGNORE_PSEUDO_BIDM;
188 
189  return IGNORE_NONE;
190 }
191 
192 
203 {
204  if(param_globals::floating_ground)
205  return;
206 
207  Stimulus* s = param_globals::stimulus;
208 
209  for(int i=0; i < param_globals::num_stim; i++) {
210  if(s[i].stimtype == Extracellular_Ground ||
211  s[i].stimtype == Extracellular_V ||
212  s[i].stimtype == Extracellular_V_OL)
213  return;
214  }
215 
216  // for now we only warn, although we should actually stop the run
217  log_msg( NULL, 4, 0,"Elliptic system is singular!\n"
218  "Either set floating_ground=1 or use an explicit ground:voltage (stimulus[X].stimtype=3)\n"
219  "Do not trust the elliptic solution of this simulation run!\n");
220 }
221 
223 {
224  printf("\n""*** GIT tag: %s\n", GIT_COMMIT_TAG);
225  printf( "*** GIT hash: %s\n", GIT_COMMIT_HASH);
226  printf( "*** GIT repo: %s\n", GIT_PATH);
227  printf( "*** dependency commits: %s\n\n", SUBREPO_COMMITS);
228 }
229 
231 {
232  // convert time steps to milliseconds
233  param_globals::dt /= 1000.;
234 
235  // check parab-solve solution method
236  if(param_globals::mass_lumping == 0 && param_globals::parab_solve==0) {
237  log_msg(NULL, 2, ECHO,
238  "Warning: explicit solve not possible without mass lumping. \n"
239  "Switching to Crank-Nicolson!\n\n");
240 
241  param_globals::parab_solve = 1;
242  }
243 
244  // check if we have to modify stimuli based on used bidomain setting
245  if(!param_globals::extracell_monodomain_stim)
246  ignore_extracellular_stim(param_globals::stimulus, param_globals::num_stim,
247  set_ignore_flags(param_globals::bidomain));
248 
249  // check nullspace if necessary
250  // if((param_globals::bidomain==BIDOMAIN && !param_globals::ellip_use_pt) ||
251  // (param_globals::bidomain==PSEUDO_BIDM && !param_globals::parab_use_pt))
252  // check_nullspace_ok();
253 
254  if(param_globals::t_sentinel > 0 && param_globals::sentinel_ID < 0 ) {
255  log_msg(0,4,0, "Warning: t_sentinel is set but no sentinel_ID has been specified; check_quiescence() behavior may not be as expected");
256  }
257 
258  if(param_globals::num_external_imp > 0 ) {
259  for(int ext_imp_i = 0; ext_imp_i < param_globals::num_external_imp; ext_imp_i++) {
260  if(param_globals::external_imp[ext_imp_i][0] != '/') {
261  log_msg(0,5,0, "external_imp[%d] error: absolute paths must be used for .so file loading (\'%s\')",
262  ext_imp_i, param_globals::external_imp[ext_imp_i]);
263  EXIT(1);
264  }
265  }
266  }
267 
268  if(param_globals::experiment == EXP_LAPLACE && param_globals::bidomain != 1) {
269  log_msg(0,4,0, "Warning: Laplace experiment mode requires bidomain = 1. Setting bidomain = 1.");
270  param_globals::bidomain = 1;
271  }
272 
273  if(param_globals::num_phys_regions == 0) {
274  log_msg(0,4,0, "Warning: No physics region defined! Please set phys_region parameters to correctly define physics.");
275 
276  if(param_globals::experiment != EXP_LAPLACE) {
277  log_msg(0,4,0, "Intra-elec and Extra-elec domains will be derived from fibers.\n");
278  param_globals::num_phys_regions = param_globals::bidomain ? 2 : 1;
279  param_globals::phys_region = (p_region*) malloc(param_globals::num_phys_regions * sizeof(p_region));
280  param_globals::phys_region[0].ptype = PHYSREG_INTRA_ELEC;
281  param_globals::phys_region[0].name = strdup("Autogenerated intracellular Electrics");
282  param_globals::phys_region[0].num_IDs = 0;
283 
284  if(param_globals::bidomain) {
285  param_globals::phys_region[1].ptype = PHYSREG_EXTRA_ELEC;
286  param_globals::phys_region[1].name = strdup("Autogenerated extracellular Electrics");
287  param_globals::phys_region[1].num_IDs = 0;
288  }
289  } else {
290  log_msg(0,4,0, "Laplace domain will be derived from fibers.\n");
291  param_globals::num_phys_regions = 1;
292  param_globals::phys_region = (p_region*) malloc(param_globals::num_phys_regions * sizeof(p_region));
293  param_globals::phys_region[0].ptype = PHYSREG_LAPLACE;
294  param_globals::phys_region[0].name = strdup("Autogenerated Laplace");
295  param_globals::phys_region[0].num_IDs = 0;
296  }
297  }
298 
299  if(param_globals::experiment == EXP_LAPLACE && !phys_defined(PHYSREG_LAPLACE)) {
300  log_msg(0,4,0, "Warning: Laplace experiment mode requires a laplace physics regions defined.");
301 
302  int idx = -1;
303  if((idx = get_phys_index(PHYSREG_EXTRA_ELEC)) > -1) {
304  log_msg(0,4,0, "Converting the defined extracellular-electrics-region to laplace-region.");
305  param_globals::phys_region[idx].ptype = PHYSREG_LAPLACE;
306  } else if ((idx = get_phys_index(PHYSREG_INTRA_ELEC)) > -1) {
307  log_msg(0,4,0, "Converting the defined intracellular-electrics-region to laplace-region.");
308  param_globals::phys_region[idx].ptype = PHYSREG_LAPLACE;
309  } else {
310  param_globals::num_phys_regions += 1;
311  param_globals::phys_region = (p_region*) realloc(param_globals::phys_region, param_globals::num_phys_regions * sizeof(p_region));
312 
313  param_globals::phys_region[param_globals::num_phys_regions - 1].ptype = PHYSREG_LAPLACE;
314  param_globals::phys_region[param_globals::num_phys_regions - 1].name = strdup("Autogenerated Laplace");
315  param_globals::phys_region[param_globals::num_phys_regions - 1].num_IDs = 0;
316  }
317  }
318 
319 #ifndef WITH_PARMETIS
320  if(param_globals::pstrat == 1) {
321  log_msg(0,3,0, "openCARP was built without Parmetis support. Swithing to KDtree.");
322  param_globals::pstrat = 2;
323  }
324 #endif
325 
326  // check if we have the legacy stimuli or the new stimuli defined by the user
327  bool legacy_stim_set = false, new_stim_set = false;
328 
329  for(int i=0; i<param_globals::num_stim; i++) {
330  Stimulus & legacy_stim = param_globals::stimulus[i];
331  Stim & new_stim = param_globals::stim[i];
332 
333  if(legacy_stim.stimtype || legacy_stim.strength)
334  legacy_stim_set = true;
335 
336  if(new_stim.crct.type || new_stim.pulse.strength)
337  new_stim_set = true;
338  }
339 
340  if(legacy_stim_set || new_stim_set) {
341  if(legacy_stim_set && new_stim_set) {
342  log_msg(0,4,0, "Warning: Legacy stimuli and default stimuli are defined. Only default stimuli will be used!");
343  }
344  else if (legacy_stim_set) {
345  log_msg(0,1,0, "Warning: Legacy stimuli defined. Please consider switching to stimulus definition \"stim[]\"!");
347  }
348  }
349  else {
350  log_msg(0,4,0, "Warning: No potential or current stimuli found!");
351  }
352 }
353 
354 void set_io_dirs(char *sim_ID, char *pp_ID, IO_t init)
355 {
356  int flg = 0, err = 0, rank = get_rank();
357 
358  char *ptr = getcwd(current_dir, 1024);
359  if (ptr == NULL) err++;
360  ptr = getcwd(input_dir, 1024);
361  if (ptr == NULL) err++;
362  //if (param_globals::experiment == 4 && post_processing_opts == MECHANIC_POSTPROCESS)
363  // { sim_ID = param_globals::ppID; param_globals::ppID = "POSTPROC_DIR"; }
364 
365  // output directory
366  if (rank == 0) {
367  if (strcmp(sim_ID, "OUTPUT_DIR")) {
368  if (mkdir(sim_ID, 0775)) { // rwxrwxr-x
369  if (errno == EEXIST ) {
370  log_msg(NULL, 2, 0, "Output directory exists: %s\n", sim_ID);
371  } else {
372  log_msg(NULL, 5, 0, "Unable to make output directory\n");
373  flg = 1;
374  }
375  }
376  } else if (mkdir(sim_ID, 0775) && errno != EEXIST) {
377  log_msg(NULL, 5, 0, "Unable to make output directory\n");
378  flg = 1;
379  }
380  }
381 
382  // terminate?
383  if(get_global(flg, MPI_SUM)) { EXIT(-1); }
384 
385  err += chdir(sim_ID);
386  ptr = getcwd(output_dir, 1024);
387  if (ptr == NULL) err++;
388 
389  // terminate?
390  if(get_global(err, MPI_SUM)) { EXIT(-1); }
391 
392  err += chdir(output_dir);
393 
394  // postprocessing directory
395  if (rank == 0 && (param_globals::experiment==EXP_POSTPROCESS)) {
396 
397  if (strcmp(param_globals::ppID, "POSTPROC_DIR")) {
398  if (mkdir(param_globals::ppID, 0775)) { // rwxrwxr-x
399  if (errno == EEXIST ) {
400  log_msg(NULL, 2, ECHO, "Postprocessing directory exists: %s\n\n", param_globals::ppID);
401  } else {
402  log_msg(NULL, 5, ECHO, "Unable to make postprocessing directory\n\n");
403  flg = 1;
404  }
405  }
406  } else if (mkdir(param_globals::ppID, 0775) && errno != EEXIST) {
407  log_msg(NULL, 5, ECHO, "Unable to make postprocessing directory\n\n");
408  flg = 1;
409  }
410 
411  }
412 
413  if(get_global(flg, MPI_SUM)) { EXIT(-1); }
414 
415  err += chdir(param_globals::ppID);
416  ptr = getcwd(postproc_dir, 1024);
417  if (ptr == NULL) err++;
418  err = chdir(output_dir);
419  if(get_global(err, MPI_SUM)) { EXIT(-1); }
420 
421  err = set_dir(init);
422  if(get_global(err, MPI_SUM)) { EXIT(-1); }
423 }
424 
425 bool setup_IO(int argc, char **argv)
426 {
427  bool io_node = false;
428  int psize = get_size(), prank = get_rank();
429 
430  if (param_globals::num_io_nodes > 0) {
431  // Can't do async IO with only one core
432  if (get_size() == 1) {
433  log_msg(NULL, 5, 0, "You cannot run with async IO on only one core.\n");
434  EXIT(EXIT_FAILURE);
435  }
436  // Can't do async IO with more IO cores than compute cores
437  if (2 * param_globals::num_io_nodes >= psize) {
438  log_msg(NULL, 5, 0, "The number of IO cores be less " "than the number of compute cores.");
439  EXIT(EXIT_FAILURE);
440  }
441 #if 0
442  if (param_globals::num_PS_nodes && param_globals::num_io_nodes > param_globals::num_PS_nodes) {
443  LOG_MSG(NULL, 5, 0,
444  "The number of IO cores (%d) should not "
445  "exceed the number of PS compute cores (%d).\n",
446  param_globals::num_io_nodes, param_globals::num_PS_nodes);
447  EXIT(-1);
448  }
449 #endif
450  // root IO node is global node 0
451  io_node = prank < param_globals::num_io_nodes;
452 
453  MPI_Comm comm;
454  MPI_Comm_split(PETSC_COMM_WORLD, io_node, get_rank(), &comm);
455  MPI_Comm_set_name(comm, io_node ? "IO" : "compute");
456 
457  PETSC_COMM_WORLD = comm; // either the compute world or IO world
458 
459  prank = get_rank();
460 
461  MPI_Intercomm_create(comm, 0, MPI_COMM_WORLD, io_node ? param_globals::num_io_nodes : 0,
463 
465  log_msg(NULL, 4, 0, "Global node %d, Comm rank %d != Intercomm rank %d\n",
466  get_rank(MPI_COMM_WORLD), get_rank(PETSC_COMM_WORLD),
468  } else
469  MPI_Comm_set_name(PETSC_COMM_WORLD, "compute");
470 
471  set_io_dirs(param_globals::simID, param_globals::ppID, OUTPUT);
472 
473  if((io_node || !param_globals::num_io_nodes) && !prank)
474  output_parameter_file("parameters.par", argc, argv);
475 
476  return io_node;
477 }
479 {
480  getcwd(current_dir, 1024);
481 }
482 
483 int set_dir(IO_t dest)
484 {
485  int err;
486 
487  if (dest==OUTPUT) err = chdir(output_dir);
488  else if (dest==POSTPROC) err = chdir(postproc_dir);
489  else if (dest==CURDIR) err = chdir(current_dir);
490  else err = chdir(input_dir);
491 
492  return err;
493 }
494 
496 {
497  // if we restart from a checkpoint, the timer_manager will be notified at a later stage
498  double start_time = 0.0;
499  user_globals::tm_manager = new timer_manager(param_globals::dt, start_time, param_globals::tend);
500 
501  double end_time = param_globals::tend;
503 
504  if(param_globals::experiment == EXP_LAPLACE) {
505  tm.initialize_singlestep_timer(tm.time, 0, iotm_console, "IO (console)", nullptr);
506  tm.initialize_singlestep_timer(tm.time, 0, iotm_state_var, "IO (state vars)", nullptr);
507  tm.initialize_singlestep_timer(tm.time, 0, iotm_spacedt, "IO (spacedt)", nullptr);
508  }
509  else {
510  tm.initialize_eq_timer(tm.time, end_time, 0, param_globals::timedt, 0, iotm_console, "IO (console)");
511  tm.initialize_eq_timer(tm.time, end_time, 0, param_globals::spacedt, 0, iotm_state_var, "IO (state vars)");
512  tm.initialize_eq_timer(tm.time, end_time, 0, param_globals::spacedt, 0, iotm_spacedt, "IO (spacedt)");
513  }
514 
515  if(param_globals::num_tsav) {
516  std::vector<double> trig(param_globals::num_tsav);
517  for(size_t i=0; i<trig.size(); i++) trig[i] = param_globals::tsav[i];
518 
519  tm.initialize_neq_timer(trig, 0, iotm_chkpt_list, "instance checkpointing");
520  }
521 
522  if(param_globals::chkpt_intv)
523  tm.initialize_eq_timer(param_globals::chkpt_start, param_globals::chkpt_stop, 0,
524  param_globals::chkpt_intv, 0, iotm_chkpt_intv, "interval checkpointing");
525 
526  if(param_globals::num_trace)
527  tm.initialize_eq_timer(tm.time, end_time, 0, param_globals::tracedt, 0, iotm_trace, "IO (node trace)");
528 }
529 
530 void get_protocol_column_widths(std::vector<int> & col_width, std::vector<int> & used_timer_ids)
531 {
532  char buff[256];
533  const short padding = 4;
534  Electrics* elec = (Electrics*) get_physics(elec_phys, false);
535 
536  do {
537  snprintf(buff, sizeof buff, "%.3lf", user_globals::tm_manager->time);
538  if(col_width[0] < int(strlen(buff)+padding))
539  col_width[0] = strlen(buff)+padding;
540 
541  snprintf(buff, sizeof buff, "%.3d", user_globals::tm_manager->d_time);
542  if(col_width[1] < int(strlen(buff)+padding))
543  col_width[1] = strlen(buff)+padding;
544 
545  int col = 2;
546  for (size_t tid = 0; tid < used_timer_ids.size(); tid++)
547  {
548  int timer_id = used_timer_ids[tid];
550 
551  if(t->d_trigger_dur && elec) {
552  // figure out value of signal linked to this timer
553  double val = 0.;
554 
555  // determine timer linked to which physics, for now we deal with electrics only
556  val = elec->timer_val(timer_id);
557 
558  snprintf(buff, sizeof buff, "%.3lf", val);
559  if(col_width[col] < int(strlen(buff)+padding))
560  col_width[col] = strlen(buff)+padding;
561  }
562  col++;
563  }
564 
565  // advance time
567  } while (!user_globals::tm_manager->elapsed());
568 
570 }
573 int plot_protocols(const char *fname)
574 {
575  int err = {0};
576  std::ofstream fh;
577  const char* smpl_endl = "\n";
578 
579  if(!get_rank()) {
580  fh.open(fname);
581 
582  // If we couldn't open the output file stream for writing
583  if (!fh) {
584  // Print an error and exit
585  log_msg(0,5,0,"Protocol file %s could not be opened for writing!\n", fname);
586  err = -1;
587  }
588  }
589 
590  // broadcast and return if err
591  if(get_global(err, MPI_SUM))
592  return err;
593 
594  // only rank 0 writes
595  if(!get_rank()) {
596 
597  // collect timer information, label, short label, unit
598  std::vector<std::string> col_labels = {"time", "tick"};
599  std::vector<std::string> col_short_labels = {"A", "B"};
600  std::vector<std::string> col_unit_labels = {"ms", "--" };
601  std::vector<int> col_width = {4, 4};
602 
603  char c_label = {'C'};
604  std::string label = {""};
605  std::string unit = {""};
606 
607  // here we store the IDs of the timers that we care about. currently this are the IO and TS timers
608  // and the electricts timers
609  std::vector<int> used_timer_ids;
610  std::vector<int> used_stim_ids;
611 
612  Electrics* elec = (Electrics*) get_physics(elec_phys, false);
613  if(elec) {
614  int sidx = 0;
615  for(const stimulus & s : elec->stimuli) {
616  if(s.ptcl.timer_id > -1) {
618  if(t) {
619  used_timer_ids.push_back(s.ptcl.timer_id);
620  used_stim_ids.push_back(sidx);
621  }
622  }
623 
624  sidx++;
625  }
626  }
627 
628  // determine longest timer label
629  int mx_llen = 0;
630  for (size_t tid = 0; tid < used_timer_ids.size(); tid++)
631  {
632  int timer_id = used_timer_ids[tid];
634 
635  int llen = strlen(t->name);
636  mx_llen = llen > mx_llen ? llen : mx_llen;
637  }
638 
639  for (size_t tid = 0; tid < used_timer_ids.size(); tid++)
640  {
641  int timer_id = used_timer_ids[tid];
643 
644  col_labels.push_back(t->name);
645  label = c_label;
646  col_short_labels.push_back(label);
647 
648  if(elec) {
649  // search physics for signals linked to timer
650  unit = elec->timer_unit(timer_id);
651  if(unit.empty()) unit = "--";
652  col_unit_labels.push_back(unit);
653  col_width.push_back(4);
654  }
655  c_label++;
656  }
657 
658  get_protocol_column_widths(col_width, used_timer_ids);
659 
660  // print header + legend first
661  fh << "# Protocol header\n#\n" << "# Legend:\n";
662  for(size_t i = 0; i<col_short_labels.size(); i++)
663  {
664  fh << "#" << std::setw(2) << col_short_labels[i] << " = " << std::setw(mx_llen) << col_labels[i];
665  fh << " [" << std::setw(10) << col_unit_labels[i] << "]";
666 
667  if(i >= 2 && used_stim_ids[i-2] > -1) {
668  stimulus & s = elec->stimuli[used_stim_ids[i-2]];
669 
670  if (is_voltage(s.phys.type)) {
671  if(s.phys.type == GND_ex)
672  fh << " ground stim" << smpl_endl;
673  else
674  fh << " applied: " << std::to_string(s.pulse.strength) << smpl_endl;
675  } else {
676  fh << smpl_endl;
677  }
678  } else {
679  fh << smpl_endl;
680  }
681  }
682 
683  // plot column short labels
684  fh << "#";
685  for(size_t i = 0; i<col_short_labels.size(); i++)
686  fh << std::setw(col_width[i] - 3) << col_short_labels[i].c_str() << std::setw(3) << " ";
687 
688  // plot column units
689  fh << smpl_endl << "#";
690  for(size_t i = 0; i<col_unit_labels.size(); i++)
691  fh << "[" << std::setw(col_width[i]-2) << col_unit_labels[i].c_str() << "]";
692 
693  // step through simulated time period
694  fh << smpl_endl << std::fixed;
695  do {
696  // time and discrete time
697  fh << std::setw(col_width[0]) << std::setprecision(3) << user_globals::tm_manager->time;
698  fh << std::setw(col_width[1]) << user_globals::tm_manager->d_time;
699 
700  // iterate over all timers
701  int col = 2;
702  for (size_t tid = 0; tid < used_timer_ids.size(); tid++)
703  {
704  int timer_id = used_timer_ids[tid];
706 
707  // type of timer: plain trigger or trigger linked to signal
708  if(!t->d_trigger_dur) {
709  int On = t->triggered ? 1 : 0;
710  fh << std::setw(col_width[col]) << On;
711  } else if(elec) {
712  // figure out value of signal linked to this timer
713  double val = 0.;
714 
715  // determine timer linked to which physics, for now we deal with electrics only
716  val = elec->timer_val(timer_id);
717 
718  fh << std::setw(col_width[col]) << std::setprecision(3) << val;
719  }
720  col++;
721  }
722 
723  fh << smpl_endl;
724 
725  // advance time
727  } while (!user_globals::tm_manager->elapsed());
728 
729  fh.close();
730 
731  // reset timer to start before actual simulation
733  }
734 
735  return err;
736 }
737 
739 {
740  const char* h1_prog = "PROG\t----- \t----\t-------\t-------|";
741  const char* h2_prog = "time\t%%comp\ttime\t ctime \t ETA |";
742  const char* h1_wc = "\tELAPS |";
743  const char* h2_wc = "\twc |";
744 
745  p.start = get_time();
746  p.last = p.start;
747 
748  log_msg(NULL, 0, 0, "%s", h1_prog );
749  log_msg(NULL, 0, NONL, "%s", h2_prog );
750  log_msg(NULL, 0, 0, "" );
751 }
752 
753 
754 void time_to_string(float time, char* str, short str_size)
755 {
756  int req_hours = ((int)(time)) / 3600;
757  int req_min = (((int)(time)) % 3600) / 60;
758  int req_sec = (((int)(time)) % 3600) % 60;
759 
760  snprintf(str, str_size, "%d:%02d:%02d", req_hours, req_min, req_sec);
761 }
762 
764 {
765 
766  float progress = 100.*(tm.time - tm.start) / (tm.end - tm.start);
767  float elapsed_time = timing(p.curr, p.start);
768  float req_time = (elapsed_time / progress) * (100.0f - progress);
769 
770  if(progress == 0.0f)
771  req_time = 0.0f;
772 
773  char elapsed_time_str[256];
774  char req_time_str[256];
775  time_to_string(elapsed_time, elapsed_time_str, 255);
776  time_to_string(req_time, req_time_str, 255);
777 
778  log_msg( NULL, 0, NONL, "%.2f\t%.1f\t%.1f\t%s\t%s",
779  tm.time,
780  progress,
781  (float)(p.curr - p.last),
782  elapsed_time_str,
783  req_time_str);
784 
785  p.last = p.curr;
786 
787  // we add an empty string for newline and flush
788  log_msg( NULL, 0, ECHO | FLUSH, "");
789 }
790 
791 void simulate()
792 {
793  // here we want to include all time dependent physics, to check if we have any of those
794  bool have_timedependent_phys = phys_defined(PHYSREG_INTRA_ELEC) || phys_defined(PHYSREG_EIKONAL);
795 
796  if(!have_timedependent_phys) {
797  log_msg(0,0,0, "\n no time-dependent physics region registered, skipping simulate loop..\n");
798  return;
799  }
800 
801  log_msg(0,0,0, "\n *** Launching simulation ***\n");
802 
803  set_dir(OUTPUT);
804 
805  if(param_globals::dump_protocol)
806  plot_protocols("protocol.trc");
807 
808  prog_stats prog;
810  init_console_output(tm, prog);
811 
812  // main loop
813  do {
814  // console output
815  if(tm.trigger(iotm_console)) {
816  // print console
817  update_console_output(tm, prog);
818  }
819 
820  // in order to be closer to carpentry we first do output and then compute the solution
821  // for the next time slice ..
822  if (tm.trigger(iotm_spacedt)) {
823  for(const auto & it : user_globals::physics_reg) {
824  it.second->output_step();
825  }
826  }
827  // compute step
828  for(const auto & it : user_globals::physics_reg) {
829  Basic_physic* p = it.second;
830  if (tm.trigger(p->timer_idx))
831  p->compute_step();
832  }
833 
834  // advance time
835  tm.update_timers();
836  } while (!tm.elapsed());
837 
838  log_msg(0,0,0, "\n\nTimings of individual physics:");
839  log_msg(0,0,0, "------------------------------\n");
840 
841  for(const auto & it : user_globals::physics_reg) {
842  Basic_physic* p = it.second;
843  p->output_timings();
844  }
845 
846 }
847 
849 {
850  if(param_globals::post_processing_opts & RECOVER_PHIE) {
851  log_msg(NULL,0,ECHO,"\nPOSTPROCESSOR: Recovering Phie ...");
852  log_msg(NULL,0,ECHO, "----------------------------------\n");
853 
854  // do postprocessing
855  int err = postproc_recover_phie();
856 
857  if(!err) {
858  log_msg(NULL,0,ECHO,"\n-----------------------------------------");
859  log_msg(NULL,0,ECHO, "POSTPROCESSOR: Successfully recoverd Phie.\n");
860  }
861  }
862 }
863 
864 Basic_physic* get_physics(physic_t p, bool error_if_missing)
865 {
866  auto it = user_globals::physics_reg.find(p);
867 
868  if(it != user_globals::physics_reg.end()) {
869  return it->second;
870  } else {
871  if(error_if_missing) {
872  log_msg(0,5,0, "%s error: required physic is not active! Usually this is due to an inconsistent experiment configuration. Aborting!", __func__);
873  EXIT(EXIT_FAILURE);
874  }
875 
876  return NULL;
877  }
878 }
879 
881 {
882  sf_vec* ret = NULL;
883 
885  ret = user_globals::datavec_reg[d];
886 
887  return ret;
888 }
889 
891 {
892  if(user_globals::datavec_reg.count(d) == 0) {
893  user_globals::datavec_reg[d] = dat;
894  }
895  else {
896  log_msg(0,5,0, "%s warning: trying to register already registered data vector.", __func__);
897  }
898 }
899 
901 {
902  std::map<mesh_t, sf_mesh> & mesh_registry = user_globals::mesh_reg;
903 
904  // This is the initial grid we read the hard-disk data into
905  mesh_registry[reference_msh] = sf_mesh();
906  // we specify the MPI communicator for the reference mesh,
907  // all derived meshes will get this comminicator automatically
908  mesh_registry[reference_msh].comm = PETSC_COMM_WORLD;
909 
910  auto register_new_mesh = [&] (mesh_t mt, int pidx) {
911  if(!mesh_registry.count(mt)) {
912  mesh_registry[mt] = sf_mesh();
913  mesh_registry[mt].name = param_globals::phys_region[pidx].name;
914  }
915  return &mesh_registry[mt];
916  };
917 
918  // based on cli parameters we determine which grids need to be defined
919  for(int i=0; i<param_globals::num_phys_regions; i++)
920  {
921  sf_mesh* curmesh = NULL;
922  // register mesh type
923  switch(param_globals::phys_region[i].ptype) {
924  case PHYSREG_EIKONAL:
925  case PHYSREG_INTRA_ELEC:
926  curmesh = register_new_mesh(intra_elec_msh, i);
927  break;
928 
929  case PHYSREG_LAPLACE:
930  case PHYSREG_EXTRA_ELEC:
931  curmesh = register_new_mesh(extra_elec_msh, i);
932  break;
933 
934  default:
935  log_msg(0,5,0, "Unsupported mesh type %d! Aborting!", param_globals::phys_region[i].ptype);
936  EXIT(EXIT_FAILURE);
937  }
938 
939  if(curmesh) {
940  // set mesh unique tags
941  for(int j=0; j<param_globals::phys_region[i].num_IDs; j++)
942  curmesh->extr_tag.insert(param_globals::phys_region[i].ID[j]);
943  }
944  }
945 }
946 
957 void retag_elements(sf_mesh & mesh, TagRegion *tagRegs, int ntr)
958 {
959  if(ntr == 0) return;
960  // checkTagRegDefs(ntr, tagRegs);
961 
963 
964  for (int i=0; i<ntr; i++) {
965  tagreg_t type = tagreg_t(tagRegs[i].type);
966  SF::vector<mesh_int_t> elem_indices;
967 
968  if (type == tagreg_list)
969  read_indices(elem_indices, tagRegs[i].elemfile, ref_eidx, mesh.comm);
970  else {
971  geom_shape shape;
972  shape.type = geom_shape::shape_t(tagRegs[i].type);
973  shape.p0.x = tagRegs[i].p0[0];
974  shape.p0.y = tagRegs[i].p0[1];
975  shape.p0.z = tagRegs[i].p0[2];
976  shape.p1.x = tagRegs[i].p1[0];
977  shape.p1.y = tagRegs[i].p1[1];
978  shape.p1.z = tagRegs[i].p1[2];
979  shape.radius = tagRegs[i].radius;
980 
981  bool nodal = false;
982  indices_from_geom_shape(elem_indices, mesh, shape, nodal);
983  }
984 
985  if(get_global((long int)elem_indices.size(), MPI_SUM, mesh.comm) == 0)
986  log_msg(0,3,0,"Tag region %d is empty", i);
987 
988  for(size_t j=0; j<elem_indices.size(); j++)
989  mesh.tag[elem_indices[j]] = tagRegs[i].tag;
990  }
991 
992  // output the vector?
993  if(strlen(param_globals::retagfile))
994  {
995  update_cwd();
996  set_dir(OUTPUT);
997 
998  int dpn = 1;
999  SF::write_data_ascii(mesh.comm, ref_eidx, mesh.tag, param_globals::retagfile, dpn);
1000 
1001  // Set dir back to what is was prior to retagfile output
1002  set_dir(CURDIR);
1003  }
1004 }
1005 
1006 size_t renormalise_fibres(SF::vector<mesh_real_t> &fib, size_t l_numelem)
1007 {
1008  size_t renormalised_count = 0;
1009 
1010  // using pragma omp without global OMP control can lead to massive compute stalls,
1011  // as all cores may be already occupied by MPI, thus they become oversubscribed. Once
1012  // there is a global OMP control in place, we can activate this parallel for again. -Aurel, 20.01.2022
1013  // #pragma omp parallel for schedule(static) reduction(+ : renormalised_count)
1014  for (size_t i = 0; i < l_numelem; i++)
1015  {
1016  const mesh_real_t f0 = fib[3*i+0], f1 = fib[3*i+1], f2 = fib[3*i+2];
1017  mesh_real_t fibre_len = sqrt(f0*f0 + f1*f1 + f2*f2);
1018 
1019  if (fibre_len && fabs(fibre_len - 1) > 1e-3) {
1020  fib[3 * i + 0] /= fibre_len;
1021  fib[3 * i + 1] /= fibre_len;
1022  fib[3 * i + 2] /= fibre_len;
1023  renormalised_count++;
1024  }
1025  }
1026 
1027  return renormalised_count;
1028 }
1029 
1031 {
1032  log_msg(0,0,0, "\n *** Processing meshes ***\n");
1033 
1034  const std::string basename = param_globals::meshname;
1035  const int verb = param_globals::output_level;
1036  std::map<mesh_t, sf_mesh> & mesh_registry = user_globals::mesh_reg;
1037  assert(mesh_registry.count(reference_msh) == 1); // There must be a reference mesh
1038 
1039  set_dir(INPUT);
1040 
1041  // we always read into the reference mesh
1042  sf_mesh & ref_mesh = mesh_registry[reference_msh];
1043  MPI_Comm comm = ref_mesh.comm;
1044 
1045  int size, rank;
1046  double t1, t2, s1, s2;
1047  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
1048 
1050  SF::vector<mesh_int_t> ptsidx;
1051 
1052  // we add pointers to the meshes that need vertex cooridnates to this list
1053  std::list< sf_mesh* > ptsread_list;
1054 
1055  // read element mesh data
1056  t1 = MPI_Wtime(); s1 = t1;
1057  if(verb) log_msg(NULL, 0, 0,"Reading reference mesh: %s.*", basename.c_str());
1058 
1059  SF::read_elements(ref_mesh, basename);
1060  SF::read_points(basename, comm, pts, ptsidx);
1061 
1062  t2 = MPI_Wtime();
1063  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1064 
1065  bool check_fibre_normality = true;
1066  if (check_fibre_normality) {
1067  t1 = MPI_Wtime();
1068 
1069  // make sure that all fibre vectors have unit length
1070  size_t l_num_fixed_fib = renormalise_fibres(ref_mesh.fib, ref_mesh.l_numelem);
1071 
1072  size_t l_num_fixed_she = 0;
1073  if (ref_mesh.she.size() > 0)
1074  l_num_fixed_she = renormalise_fibres(ref_mesh.she, ref_mesh.l_numelem);
1075 
1076  unsigned long fixed[2] = {(unsigned long) l_num_fixed_fib, (unsigned long) l_num_fixed_she};
1077  MPI_Allreduce(MPI_IN_PLACE, fixed, 2, MPI_UNSIGNED_LONG, MPI_SUM, comm);
1078 
1079  if (fixed[0] + fixed[1] > 0)
1080  log_msg(NULL, 0, 0, "Renormalised %ld longitudinal and %ld sheet-transverse fibre vectors.", fixed[0], fixed[1]);
1081 
1082  t2 = MPI_Wtime();
1083  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1084  }
1085 
1086  if(param_globals::numtagreg > 0) {
1087  log_msg(0, 0, 0, "Re-tagging reference mesh");
1088 
1089  // the retagging requires vertex coordinates, as such we need to read them into
1090  // the reference mesh
1091  ptsread_list.push_back(&ref_mesh);
1092  SF::insert_points(pts, ptsidx, ptsread_list);
1093 
1094  retag_elements(ref_mesh, param_globals::tagreg, param_globals::numtagreg);
1095 
1096  // we clear the list of meshet to receive vertices
1097  ptsread_list.clear();
1098  }
1099 
1100  if(verb) log_msg(NULL, 0, 0, "Processing submeshes");
1101 
1102  for(auto it = mesh_registry.begin(); it != mesh_registry.end(); ++it) {
1103  mesh_t grid_type = it->first;
1104  sf_mesh & submesh = it->second;
1105 
1106  if(grid_type != reference_msh) {
1107  if(verb > 1) log_msg(NULL, 0, 0, "\nSubmesh name: %s", submesh.name.c_str());
1108  t1 = MPI_Wtime();
1109 
1110  if(submesh.extr_tag.size())
1111  extract_tagbased(ref_mesh, submesh);
1112  else {
1113  // all submeshes should be defined on sets of tags, for backwards compatibility
1114  // we do a fiber based intra_elec_msh extraction if no tags are provided. Also, we
1115  // could do special treatments of any other physics type here. It would defeat
1116  // the purpose of the tag-based design, though. -Aurel
1117  switch(grid_type) {
1118  case intra_elec_msh: extract_myocardium(ref_mesh, submesh); break;
1119  default: extract_tagbased(ref_mesh, submesh); break;
1120  }
1121  }
1122  t2 = MPI_Wtime();
1123  if(verb > 1) log_msg(NULL, 0, 0, "Extraction done in %f sec.", float(t2 - t1));
1124 
1125  ptsread_list.push_back(&submesh);
1126  }
1127  }
1128 
1129  // KDtree partitioning requires the coordinates to be present in the mesh data
1130  if(param_globals::pstrat == 2)
1131  SF::insert_points(pts, ptsidx, ptsread_list);
1132 
1133  for(auto it = mesh_registry.begin(); it != mesh_registry.end(); ++it)
1134  {
1135  mesh_t grid_type = it->first;
1136  sf_mesh & submesh = it->second;
1137  if(grid_type != reference_msh) {
1138  if(verb > 2) log_msg(NULL, 0, 0, "\nSubmesh name: %s", submesh.name.c_str());
1140 
1141  // generate partitioning vector
1142  t1 = MPI_Wtime();
1143  switch(param_globals::pstrat) {
1144  case 0:
1145  if(verb > 2) log_msg(NULL, 0, 0, "Using linear partitioning ..");
1146  break;
1147 
1148 #ifdef WITH_PARMETIS
1149  case 1:
1150  {
1151  if(verb > 2) log_msg(NULL, 0, 0, "Using Parmetis partitioner ..");
1152  SF::parmetis_partitioner<mesh_int_t, mesh_real_t> partitioner(param_globals::pstrat_imbalance, 2);
1153  partitioner(submesh, part);
1154  break;
1155  }
1156 #endif
1157  default:
1158  case 2: {
1159  if(verb > 2) log_msg(NULL, 0, 0, "Using KDtree partitioner ..");
1161  partitioner(submesh, part);
1162  break;
1163  }
1164  }
1165  t2 = MPI_Wtime();
1166  if(verb > 2) log_msg(NULL, 0, 0, "Partitioning done in %f sec.", float(t2 - t1));
1167 
1168  if(param_globals::pstrat > 0) {
1169  if(param_globals::gridout_p) {
1170  std::string out_name = get_basename(param_globals::meshname);
1171  if(grid_type == intra_elec_msh) out_name += "_i.part.dat";
1172  else if(grid_type == extra_elec_msh) out_name += "_e.part.dat";
1173 
1174  set_dir(OUTPUT);
1175  log_msg(0,0,0, "Writing \"%s\" partitioning to: %s", submesh.name.c_str(), out_name.c_str());
1176  write_data_ascii(submesh.comm, submesh.get_numbering(SF::NBR_ELEM_REF), part, out_name);
1177  }
1178 
1179  t1 = MPI_Wtime();
1180  SF::redistribute_elements(submesh, part);
1181  t2 = MPI_Wtime();
1182  if(verb > 2) log_msg(NULL, 0, 0, "Redistributing done in %f sec.", float(t2 - t1));
1183  }
1184 
1185  t1 = MPI_Wtime();
1187  sm_numbering(submesh);
1188  t2 = MPI_Wtime();
1189  if(verb > 2) log_msg(NULL, 0, 0, "Canonical numbering done in %f sec.", float(t2 - t1));
1190 
1191  t1 = MPI_Wtime();
1192  submesh.generate_par_layout();
1193  SF::petsc_numbering<mesh_int_t, mesh_real_t> p_numbering(submesh.pl);
1194  p_numbering(submesh);
1195  t2 = MPI_Wtime();
1196  if(verb > 2) log_msg(NULL, 0, 0, "PETSc numbering done in %f sec.", float(t2 - t1));
1197  if(verb > 2) print_DD_info(submesh);
1198  }
1199  }
1200 
1201  SF::insert_points(pts, ptsidx, ptsread_list);
1202  ref_mesh.clear_data();
1203 
1204  s2 = MPI_Wtime();
1205  if(verb) log_msg(NULL, 0, 0, "All done in %f sec.", float(s2 - s1));
1206 }
1207 
1209 {
1210  bool write_intra_elec = mesh_is_registered(intra_elec_msh) && param_globals::gridout_i;
1211  bool write_extra_elec = mesh_is_registered(extra_elec_msh) && param_globals::gridout_e;
1212 
1213  set_dir(OUTPUT);
1214  std::string output_base = get_basename(param_globals::meshname);
1215 
1216  if(write_intra_elec) {
1217  sf_mesh & mesh = get_mesh(intra_elec_msh);
1218 
1219  if(param_globals::gridout_i & 1) {
1220  sf_mesh surfmesh;
1221  if(param_globals::output_level > 1)
1222  log_msg(0,0,0, "Computing \"%s\" surface ..", mesh.name.c_str());
1223  compute_surface_mesh(mesh, SF::NBR_SUBMESH, surfmesh);
1224 
1225  std::string output_file = output_base + "_i.surf";
1226  log_msg(0,0,0, "Writing \"%s\" surface: %s", mesh.name.c_str(), output_file.c_str());
1227  write_surface(surfmesh, output_file);
1228  }
1229  if(param_globals::gridout_i & 2) {
1230  bool write_binary = false;
1231 
1232  std::string output_file = output_base + "_i";
1233  log_msg(0,0,0, "Writing \"%s\" mesh: %s", mesh.name.c_str(), output_file.c_str());
1234  write_mesh_parallel(mesh, write_binary, output_file.c_str());
1235  }
1236  }
1237 
1238  if(write_extra_elec) {
1239  sf_mesh & mesh = get_mesh(extra_elec_msh);
1240 
1241  if(param_globals::gridout_e & 1) {
1242  sf_mesh surfmesh;
1243  if(param_globals::output_level > 1)
1244  log_msg(0,0,0, "Computing \"%s\" surface ..", mesh.name.c_str());
1245 
1246  compute_surface_mesh(mesh, SF::NBR_SUBMESH, surfmesh);
1247  std::string output_file = output_base + "_e.surf";
1248  log_msg(0,0,0, "Writing \"%s\" surface: %s", mesh.name.c_str(), output_file.c_str());
1249  write_surface(surfmesh, output_file);
1250  }
1251  if(param_globals::gridout_e & 2) {
1252  bool write_binary = false;
1253  std::string output_file = output_base + "_e";
1254  log_msg(0,0,0, "Writing \"%s\" mesh: %s", mesh.name.c_str(), output_file.c_str());
1255  write_mesh_parallel(mesh, write_binary, output_file.c_str());
1256  }
1257  }
1258 }
1259 
1260 [[noreturn]] void cleanup_and_exit()
1261 {
1262  destroy_physics();
1264 
1265  param_clean();
1266  PetscFinalize();
1267 
1268  // close petsc error FD
1271 
1272  exit(EXIT_SUCCESS);
1273 }
1274 
1275 char* get_file_dir(const char* file)
1276 {
1277  char* filecopy = dupstr(file);
1278  char* dir = dupstr(dirname(filecopy));
1279 
1280  free(filecopy);
1281  return dir;
1282 }
1283 
1285 {
1286  int rank = get_rank();
1287  set_dir(OUTPUT);
1288 
1289  if(rank == 0) {
1290  // we open a error log file handle and set it as petsc stderr
1291  user_globals::petsc_error_fd = fopen("petsc_err_log.txt", "w");
1292  PETSC_STDERR = user_globals::petsc_error_fd;
1293  }
1294  else {
1295  PetscErrorPrintf = PetscErrorPrintfNone;
1296  }
1297 }
1298 
1300 {
1301  const sf_mesh & mesh = get_mesh(id);
1303  short mindim = 3;
1304 
1305  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
1306  view.set_elem(eidx);
1307  short cdim = view.dimension();
1308  if(mindim < cdim) mindim = cdim;
1309  }
1310 
1311  mindim = get_global(mindim, MPI_MIN, mesh.comm);
1312 
1313  return mindim;
1314 }
1315 
1317  const mesh_t inp_meshid,
1318  const int dpn,
1319  const char* name,
1320  const char* units,
1321  const SF::vector<mesh_int_t>* idx,
1322  bool elem_data)
1323 {
1324  sync_io_item IO;
1325 
1326  IO.data = inp_data;
1327  IO.elem_flag = elem_data;
1328  IO.restr_idx = idx;
1329 
1330  IGBheader regigb;
1332  const int num_io = tm.timers[iotm_spacedt]->numIOs;
1333  int err = 0;
1334 
1335  int gsize = inp_data->gsize();
1336 
1337  // if we are restricting, we have to compute the restricted global size
1338  if(idx != NULL)
1339  gsize = get_global(int(idx->size()), MPI_SUM);
1340 
1341  regigb.x(gsize / dpn);
1342  regigb.dim_x(regigb.x()-1);
1343  regigb.inc_x(1);
1344 
1345  regigb.y(1); regigb.z(1);
1346  regigb.t(num_io);
1347  regigb.dim_t(tm.end-tm.start);
1348 
1349  switch(dpn) {
1350  default:
1351  case 1: regigb.type(IGB_FLOAT); break;
1352  case 3: regigb.type(IGB_VEC3_f); break;
1353  case 4: regigb.type(IGB_VEC4_f); break;
1354  case 9: regigb.type(IGB_VEC9_f); break;
1355  }
1356 
1357  regigb.unites_x("um"); regigb.unites_y("um"); regigb.unites_z("um");
1358  regigb.unites_t("ms");
1359  regigb.unites(units);
1360 
1361  regigb.inc_t(param_globals::spacedt);
1362 
1363  if(get_rank() == 0) {
1364  FILE_SPEC file = f_open(name, "w");
1365  if(file != NULL) {
1366  regigb.fileptr(file->fd);
1367  regigb.write();
1368  delete file;
1369  }
1370  else err++;
1371  }
1372 
1373  err = get_global(err, MPI_SUM);
1374  if(err) {
1375  log_msg(0,5,0, "%s error: Could not set up data output! Aborting!", __func__);
1376  EXIT(1);
1377  }
1378 
1379  IO.igb = regigb;
1380 
1381  SF::mixed_tuple<mesh_t, int> mesh_spec = {inp_meshid, dpn};
1382  IO.spec = mesh_spec;
1383 
1384  if(elem_data) {
1385  if(buffmap_elem.find(mesh_spec) == buffmap_elem.end()) {
1386  sf_vec *inp_copy; SF::init_vector(&inp_copy, inp_data);
1387  buffmap_elem[mesh_spec] = inp_copy;
1388  }
1389  } else {
1390  if(buffmap.find(mesh_spec) == buffmap.end()) {
1391  sf_vec *inp_copy; SF::init_vector(&inp_copy, inp_data);
1392  buffmap[mesh_spec] = inp_copy;
1393  }
1394  }
1395 
1396  this->sync_IOs.push_back(IO);
1397 }
1398 
1399 void igb_output_manager::register_output_async(sf_vec* inp_data,
1400  const mesh_t inp_meshid,
1401  const int dpn,
1402  const char* name,
1403  const char* units,
1404  const SF::vector<mesh_int_t>* idx,
1405  bool elem_data)
1406 {
1407  sf_mesh & mesh = get_mesh(inp_meshid);
1408  SF::vector<mesh_int_t> ioidx;
1409  int rank = get_rank();
1410 
1411  async_io_item IO;
1412  IO.data = inp_data;
1413  IO.restr_idx = idx;
1414 
1415  if(elem_data) {
1417  ioidx.resize(mesh.l_numelem);
1418  for(size_t i=0; i<mesh.l_numelem; i++)
1419  ioidx[i] = nbr[i];
1420  } else {
1421  const SF::vector<mesh_int_t> & alg_nod = mesh.pl.algebraic_nodes();
1423 
1424  if(idx == NULL) {
1425  ioidx.resize(alg_nod.size());
1426 
1427  for(size_t i=0; i<alg_nod.size(); i++)
1428  ioidx[i] = nbr[alg_nod[i]];
1429  } else {
1430  ioidx.resize(idx->size());
1431  IO.restr_petsc_idx.resize(idx->size());
1432 
1433  for(size_t i=0; i<idx->size(); i++) {
1434  mesh_int_t loc_nodal = (*idx)[i];
1435  ioidx[i] = nbr[loc_nodal];
1436  IO.restr_petsc_idx[i] = local_nodal_to_local_petsc(mesh, rank, loc_nodal);
1437  }
1438  }
1439  }
1440 
1441  int id = async::COMPUTE_register_output(ioidx, dpn, name, units);
1442  IO.IO_id = id;
1443 
1444  this->async_IOs.push_back(IO);
1445 }
1446 
1448  const mesh_t inp_meshid,
1449  const int dpn,
1450  const char* name,
1451  const char* units,
1452  const SF::vector<mesh_int_t>* idx,
1453  bool elem_data)
1454 {
1455  if(param_globals::num_io_nodes == 0)
1456  register_output_sync(inp_data, inp_meshid, dpn, name, units, idx, elem_data);
1457  else
1458  register_output_async(inp_data, inp_meshid, dpn, name, units, idx, elem_data);
1459 }
1460 
1461 sf_vec* igb_output_manager::fill_output_buffer(const sync_io_item & it)
1462 {
1463  const SF::mixed_tuple<mesh_t, int> & cspec = it.spec;
1464  sf_vec* data_vec = it.data;
1465 
1466  bool have_perm = it.elem_flag ? have_permutation(cspec.v1, ELEM_PETSC_TO_CANONICAL, cspec.v2):
1467  have_permutation(cspec.v1, PETSC_TO_CANONICAL, cspec.v2);
1468 
1469  if(have_perm) {
1470  sf_vec* perm_vec = it.elem_flag ? this->buffmap_elem[cspec] : this->buffmap[cspec];
1472  get_permutation(cspec.v1, PETSC_TO_CANONICAL, cspec.v2);
1473  sc->forward(*data_vec, *perm_vec);
1474  return perm_vec;
1475  } else {
1476  return data_vec;
1477  }
1478 }
1479 
1481 {
1482  SF::vector<float> restr_buff;
1483  int rank = get_rank();
1484  // loop over registered datasets and root-write one by one
1485  //
1486  for (sync_io_item & it : sync_IOs) {
1487  // write to associated file descriptor
1488  FILE* fd = static_cast<FILE*>(it.igb.fileptr());
1489 
1490  // fill the output buffer
1491  sf_vec* buff = fill_output_buffer(it);
1492 
1493  if(it.restr_idx == NULL) {
1494  buff->write_binary<float>(fd);
1495  } else {
1496  const SF::vector<mesh_int_t> & idx = *it.restr_idx;
1497  SF_real* p = buff->ptr();
1498 
1499  restr_buff.resize(idx.size()); restr_buff.resize(0);
1500 
1501  for(mesh_int_t ii : idx)
1502  restr_buff.push_back(p[ii]);
1503 
1504  root_write(fd, restr_buff, PETSC_COMM_WORLD);
1505  buff->release_ptr(p);
1506  }
1507  }
1508 
1509  // do all A-synchronous output:
1510  // loop over IDs received from the IO nodes and trigger async output
1511  //
1512  for (async_io_item & it : async_IOs) {
1513  SF_real* p = it.data->ptr();
1514  int ls = it.data->lsize();
1515  int id = it.IO_id;
1516 
1517  if(it.restr_idx == NULL)
1518  async::COMPUTE_do_output(p, ls, id);
1519  else {
1520  async::COMPUTE_do_output(p, it.restr_petsc_idx, id);
1521  }
1522 
1523  it.data->release_ptr(p);
1524  }
1525 }
1526 
1528 {
1529  if(get_rank() == 0) {
1530  // loop over registered datasets and close fd
1531  for(sync_io_item & it : sync_IOs) {
1532  FILE* fd = static_cast<FILE*>(it.igb.fileptr());
1533  fclose(fd);
1534  }
1535  }
1536 
1537  for(auto it = buffmap.begin(); it != buffmap.end(); ++it)
1538  delete it->second;
1539 
1540  for(auto it = buffmap_elem.begin(); it != buffmap_elem.end(); ++it)
1541  delete it->second;
1542 
1543  // we resize the arrays and clear the maps so that we are safe when calling
1544  // close_files_and_cleanup multiple times.
1545  sync_IOs.resize(0);
1546  async_IOs.resize(0);
1547  buffmap.clear();
1548  buffmap_elem.clear();
1549 }
1550 
1552 {
1553  for(sync_io_item & it : sync_IOs) {
1554  if(it.data == vec)
1555  return &it.igb;
1556  }
1557 
1558  return NULL;
1559 }
1560 
1561 void output_parameter_file(const char *fname, int argc, char **argv)
1562 {
1563  const int max_line_len = 128;
1564  const char* file_sep = "#=======================================================";
1565 
1566  // make sure only root executes this function
1567  if(get_rank() != 0)
1568  return;
1569 
1570  FILE_SPEC out = f_open(fname, "w");
1571  fprintf(out->fd, "# CARP GIT commit hash: %s\n", GIT_COMMIT_HASH);
1572  fprintf(out->fd, "# dependency hashes: %s\n", SUBREPO_COMMITS);
1573  fprintf(out->fd, "\n");
1574 
1575  // output the command line
1576  char line[8196] = "# ";
1577 
1578  for (int j=0; j<argc; j++) {
1579  strcat(line, argv[j]);
1580  if(strlen(line) > max_line_len) {
1581  fprintf(out->fd, "%s\n", line);
1582  strcpy(line, "# ");
1583  } else
1584  strcat(line, " ");
1585  }
1586 
1587  fprintf(out->fd, "%s\n\n", line);
1588  set_dir(INPUT);
1589 
1590  // convert command line to a par file
1591  for (int i=1; i<argc; i++) {
1592  std::string argument(argv[i]);
1593  if (argument == "+F" || argument.find("_options_file")!= std::string::npos) {
1594 
1595  std::string init = "";
1596  if (argument.find("_options_file")!= std::string::npos) {
1597  fprintf(out->fd, "%s = %s\n", argument.substr(1).c_str(), argv[i+1]);
1598  init = "#";
1599  }
1600  fprintf(out->fd, "%s>>\n", file_sep);
1601  // import par files
1602  i++;
1603  fprintf(out->fd, "## %s ##\n", argv[i]);
1604  FILE *in = fopen(argv[i], "r");
1605  while (fgets(line, 8196, in))
1606  fprintf(out->fd, "%s%s", init.c_str(), line);
1607  fclose(in);
1608  fprintf(out->fd, "\n##END of %s\n", argv[i]);
1609  fprintf(out->fd, "%s<<\n\n", file_sep);
1610  }
1611  else if(argv[i][0] == '-')
1612  {
1613  bool prelast = (i==argc-1);
1614  bool paramFollows = !prelast && ((argv[i+1][0] != '-') ||
1615  ((argv[i+1][1] >= '0') && (argv[i+1][1] <= '9')));
1616 
1617  // strip leading hyphens from command line opts
1618  // assume options do not start with numbers
1619  if(paramFollows) {
1620  // nonflag option
1621  char *optcpy = strdup(argv[i+1]);
1622  char *front = optcpy;
1623  // strip {} if present for arrays of values
1624  while(*front==' ' && *front) front++;
1625  if(*front=='{') {
1626  while(*++front == ' ');
1627  char *back = optcpy+strlen(optcpy)-1;
1628  while(*back==' ' && back>front) back--;
1629  if(*back == '}')
1630  *back = '\0';
1631  }
1632  if (strstr(front, "=") != nullptr) // if "=" is find then we need ""
1633  fprintf(out->fd, "%-40s= \"%s\"\n", argv[i]+1, front);
1634  else
1635  fprintf(out->fd, "%-40s= %s\n", argv[i]+1, front);
1636  free(optcpy);
1637  i++;
1638  }
1639  else // a flag was specified
1640  fprintf(out->fd, "%-40s= 1\n", argv[i]);
1641  }
1642  }
1643  f_close(out);
1644 }
1645 
1646 void savequit()
1647 {
1648  if(!get_physics(elec_phys)) return;
1649 
1650  set_dir(OUTPUT);
1651 
1652  double time = user_globals::tm_manager->time;
1653  char save_fn[512];
1654 
1655  snprintf(save_fn, sizeof save_fn, "exit.save.%.3f.roe", time);
1656  log_msg(NULL, 0, 0, "savequit called at time %g\n", time);
1657 
1659  elec->ion.miif->dump_state(save_fn, time, intra_elec_msh, false, GIT_COMMIT_COUNT);
1660 
1661  cleanup_and_exit();
1662 }
1663 
1664 } // namespace opencarp
1665 
void generate_par_layout()
Set up the parallel layout.
Definition: SF_container.h:526
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.
void register_physics()
Register physics to the physics registry.
Definition: sim_utils.cc:89
void update_console_output(const timer_manager &tm, prog_stats &p)
Definition: sim_utils.cc:763
bool elem_flag
igb header we use for output
Definition: sim_utils.h:276
#define IGNORE_NONE
Definition: stimulate.h:54
std::map< physic_t, Basic_physic * > physics_reg
the physics
Definition: main.cc:50
IO_t
The different output (directory) types.
Definition: sim_utils.h:53
#define STM_IGNORE_MONODOMAIN
Definition: stimulate.h:60
virtual void release_ptr(S *&p)=0
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:452
class to store shape definitions
Definition: basics.h:381
void output_parameter_file(const char *fname, int argc, char **argv)
Definition: sim_utils.cc:1561
#define FLUSH
Definition: basics.h:311
stim_pulse pulse
stimulus wave form
Definition: stimulate.h:165
void set_io_dirs(char *sim_ID, char *pp_ID, IO_t init)
Definition: sim_utils.cc:354
int get_phys_index(int physreg)
get index in param_globals::phys_region array for a certain phys region
char * dupstr(const char *old_str)
Definition: basics.cc:44
void dim_t(float a)
Definition: IGBheader.h:333
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
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:392
const SF::vector< mesh_int_t > * restr_idx
when using asyncIO, here we store the different IDs associated to the vectors we output ...
Definition: sim_utils.h:283
#define IGB_FLOAT
Definition: IGBheader.h:60
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.
float mesh_real_t
Definition: SF_container.h:47
vector< T > tag
element tag
Definition: SF_container.h:405
tagreg_t
tag regions types. must be in line with carp.prm
Definition: sim_utils.h:56
PETSc numbering of nodes.
Definition: SF_container.h:191
#define PETSC_TO_CANONICAL
Permute algebraic data from PETSC to canonical ordering.
Definition: sf_interface.h:74
#define IGB_VEC9_f
Definition: IGBheader.h:76
void unites_t(const char *a)
Definition: IGBheader.h:354
Comfort class. Provides getter functions to access the mesh member variables more comfortably...
Definition: SF_fem_utils.h:704
#define BIDOMAIN
Definition: sim_utils.h:144
void time_to_string(float time, char *str, short str_size)
Definition: sim_utils.cc:754
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:417
std::map< mesh_t, sf_mesh > mesh_reg
Registry for the different meshes used in a multi-physics simulation.
Definition: main.cc:46
vector< S > fib
fiber direction
Definition: SF_container.h:408
void reset_timers()
Reset time in timer_manager and then reset registered timers.
Definition: timer_utils.h:115
void free_scatterings()
Free the registered scatterings.
bool is_voltage(stim_t type)
uses current for stimulation
Definition: stimulate.cc:65
void simulate()
Main simulate loop.
Definition: sim_utils.cc:791
void post_process()
do postprocessing
Definition: sim_utils.cc:848
void initialize_eq_timer(double istart, double iend, int ntrig, double iintv, double idur, int ID, const char *iname, const char *poolname=nullptr)
Definition: timer_utils.cc:48
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
Definition: electrics.cc:717
size_t renormalise_fibres(SF::vector< mesh_real_t > &fib, size_t l_numelem)
Definition: sim_utils.cc:1006
virtual T lsize() const =0
bool triggered
flag indicating trigger at current time step
Definition: timer_utils.h:54
void write_mesh_parallel(const meshdata< T, S > &mesh, bool binary, std::string basename)
Write a parallel mesh to harddisk without gathering it on one rank.
#define IGB_VEC3_f
Definition: IGBheader.h:71
timer_manager * tm_manager
a manager for the various physics timers
Definition: main.cc:52
std::map< SF::mixed_tuple< mesh_t, int >, sf_vec * > buffmap
map data spec -> PETSc vector buffer
Definition: sim_utils.h:308
SF::vector< mesh_int_t > restr_petsc_idx
pointer to index vector with nodal indices we restrict to.
Definition: sim_utils.h:284
double last
last output wallclock time
Definition: sim_utils.h:268
void setup_petsc_err_log()
set up error logs for PETSc, so that it doesnt print errors to stderr.
Definition: sim_utils.cc:1284
void write_data_ascii(const MPI_Comm comm, const vector< T > &idx, const vector< S > &data, std::string file, short dpn=1)
void indices_from_geom_shape(SF::vector< mesh_int_t > &idx, const sf_mesh &mesh, const geom_shape shape, const bool nodal)
Populate vertex data with the vertices inside a defined box shape.
Definition: fem_utils.cc:169
const char * name
The name of the physic, each physic should have one.
Definition: physics_types.h:62
void initialize_physics()
Initialize all physics in the registry.
Definition: sim_utils.cc:103
void savequit()
save state and quit simulator
Definition: sim_utils.cc:1646
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
#define ECHO
Definition: basics.h:308
T * data()
Pointer to the vector&#39;s start.
Definition: SF_vector.h:91
void read_elements(meshdata< T, S > &mesh, std::string basename)
Read the element data (elements and fibers) of a CARP mesh.
Definition: SF_mesh_io.h:546
IGBheader igb
pointer to index vector used for restricting output.
Definition: sim_utils.h:275
mesh_t
The enum identifying the different meshes we might want to load.
Definition: sf_interface.h:58
T & push_back(T val)
Definition: SF_vector.h:283
bool trigger(int ID) const
Definition: timer_utils.h:166
physic_t
Identifier for the different physics we want to set up.
Definition: physics_types.h:51
double end
final time
Definition: timer_utils.h:81
void close_files_and_cleanup()
close file descriptors
Definition: sim_utils.cc:1527
void check_and_convert_params()
Here we want to put all parameter checks, conversions and modifications that have been littered throu...
Definition: sim_utils.cc:230
void init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:99
#define NO_EXTRA_V
Definition: stimulate.h:56
void write_surface(const meshdata< T, S > &surfmesh, std::string surffile)
Definition: SF_mesh_io.h:732
#define IsExtraV(A)
Definition: stimulate.h:50
#define ELEM_PETSC_TO_CANONICAL
Permute algebraic element data from PETSC to canonical ordering.
Definition: sf_interface.h:76
char * get_file_dir(const char *file)
Definition: sim_utils.cc:1275
void inc_t(float a)
Definition: IGBheader.h:321
void insert_points(const vector< S > &pts, const vector< T > &ptsidx, std::list< meshdata< T, S > *> &meshlist)
Insert the points from the read-in buffers into a list of distributed meshes.
Definition: SF_mesh_io.h:915
void COMPUTE_do_output(SF_real *dat, const int lsize, const int IO_id)
Definition: async_io.cc:396
for display execution progress and statistical data of electrical solve
Definition: sim_utils.h:265
bool setup_IO(int argc, char **argv)
Definition: sim_utils.cc:425
virtual S * ptr()=0
bool using_legacy_stimuli
flag storing whether legacy stimuli are used
Definition: main.cc:58
short dimension() const
Definition: SF_fem_utils.h:920
void dump_state(char *, float, opencarp::mesh_t gid, bool, unsigned int)
#define RECOVER_PHIE
Definition: sim_utils.h:149
#define NO_EXTRA_GND
Definition: stimulate.h:55
#define EXP_LAPLACE
Definition: sim_utils.h:159
FILE * petsc_error_fd
file descriptor for petsc error output
Definition: main.cc:56
Base class for tracking progress.
Definition: progress.hpp:38
int mesh_int_t
Definition: SF_container.h:46
void initialize_singlestep_timer(double tg, double idur, int ID, const char *iname, const char *poolname=nullptr)
Definition: timer_utils.h:156
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:192
void unites(const char *a)
Definition: IGBheader.h:357
double start
output start wallclock time
Definition: sim_utils.h:267
void extract_tagbased(const meshdata< T, S > &mesh, meshdata< T, S > &submesh)
Extract a submesh based on element tags.
void fileptr(gzFile f)
Definition: IGBheader.cc:336
Submesh nodal numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:190
#define PHYSREG_INTRA_ELEC
Definition: sf_interface.h:79
virtual void destroy()=0
void unites_z(const char *a)
Definition: IGBheader.h:351
void init_console_output(const timer_manager &tm, prog_stats &p)
Definition: sim_utils.cc:738
double time
current time
Definition: timer_utils.h:76
Async IO functions.
size_t size() const
Definition: hashmap.hpp:1066
void cleanup_and_exit()
Definition: sim_utils.cc:1260
int COMPUTE_register_output(const SF::vector< mesh_int_t > &idx, const int dpn, const char *name, const char *units)
Definition: async_io.cc:104
The abstract physics interface we can use to trigger all physics.
Definition: physics_types.h:59
Top-level header of FEM module.
#define Intracellular_I
Definition: stimulate.h:41
#define PHYSREG_EXTRA_ELEC
Definition: sf_interface.h:80
void forward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Forward scattering.
void parse_mesh_types()
Parse the phys_type CLI parameters and set up (empty) SF::meshdata meshes.
Definition: sim_utils.cc:900
void parse_params_cpy(int argc, char **argv)
Initialize input parameters on a copy of the real command line parameters.
Definition: sim_utils.cc:46
int set_ignore_flags(int mode)
Definition: sim_utils.cc:180
double curr
current output wallclock time
Definition: sim_utils.h:269
void unites_x(const char *a)
Definition: IGBheader.h:345
IGBheader * get_igb_header(const sf_vec *vec)
Get the pointer to the igb header for a vector that was registered for output.
Definition: sim_utils.cc:1551
void destroy_physics()
Destroy all physics in the registry.
Definition: sim_utils.cc:127
#define Extracellular_I
Definition: stimulate.h:38
int timer_idx
the timer index received from the timer manager
Definition: physics_types.h:66
void write_data()
write registered data to disk
Definition: sim_utils.cc:1480
long d_time
current time instance index
Definition: timer_utils.h:77
const SF::vector< mesh_int_t > * restr_idx
pointer to data registered for output
Definition: sim_utils.h:274
std::map< SF::mixed_tuple< mesh_t, int >, sf_vec * > buffmap_elem
Definition: sim_utils.h:309
void dim_x(float a)
Definition: IGBheader.h:324
void basic_timer_setup()
Here we set up the timers that we always want to have, independent of physics.
Definition: sim_utils.cc:495
void clear_data()
Clear the mesh data from memory.
Definition: SF_container.h:540
sf_mesh & get_mesh(const mesh_t gt)
Get a mesh by specifying the gridID.
Definition: sf_interface.cc:33
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:1316
int plot_protocols(const char *fname)
plot simulation protocols (I/O timers, stimuli, boundary conditions, etc)
Definition: sim_utils.cc:573
Submesh element numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:193
void print_DD_info(const meshdata< T, S > &mesh)
Print some basic information on the domain decomposition of a mesh.
size_t write_binary(FILE *fd)
Write a vector to HD in binary. File descriptor is already set up.
#define NO_EXTRA_I
Definition: stimulate.h:57
void get_protocol_column_widths(std::vector< int > &col_width, std::vector< int > &used_timer_ids)
Definition: sim_utils.cc:530
void setup_meshes()
Read in the reference mesh and use its data to populate all meshes registered in the mesh registry...
Definition: sim_utils.cc:1030
#define PHYSREG_EIKONAL
Definition: sf_interface.h:81
hashmap::unordered_set< int > extr_tag
the element tags based on which the mesh has been extracted
Definition: SF_container.h:412
void inc_x(float a)
Definition: IGBheader.h:312
long d_trigger_dur
discrete duration
Definition: timer_utils.h:58
void show_build_info()
show the build info, exit if -buildinfo was provided. This code runs before MPI_Init().
Definition: sim_utils.cc:222
void insert(InputIterator first, InputIterator last)
Definition: hashmap.hpp:962
#define PSEUDO_BIDM
Definition: sim_utils.h:145
void output_meshes()
Definition: sim_utils.cc:1208
#define STM_IGNORE_BIDOMAIN
Definition: stimulate.h:59
void register_data(sf_vec *dat, datavec_t d)
Register a data vector in the global registry.
Definition: sim_utils.cc:890
void get_time(double &tm)
Definition: basics.h:436
Basic_physic * get_physics(physic_t p, bool error_if_missing)
Convinience function to get a physics.
Definition: sim_utils.cc:864
void set_elem(size_t eidx)
Set the view to a new element.
Definition: SF_fem_utils.h:732
void extract_myocardium(const meshdata< T, S > &mesh, meshdata< T, S > &submesh)
Extract the myocardium submesh.
#define EXP_POSTPROCESS
Definition: sim_utils.h:161
void retag_elements(sf_mesh &mesh, TagRegion *tagRegs, int ntr)
Definition: sim_utils.cc:957
int timer_id
timer for stimulus
Definition: stimulate.h:121
#define PHYSREG_LAPLACE
Definition: sf_interface.h:84
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
Functor class generating a numbering optimized for PETSc.
Definition: SF_numbering.h:230
SF::mixed_tuple< mesh_t, int > spec
flag whether the data is elements-wise
Definition: sim_utils.h:277
MPI_Comm IO_Intercomm
Communicator between IO and compute worlds.
Definition: main.cc:60
#define Ignore_Stim
Definition: stimulate.h:48
V timing(V &t2, const V &t1)
Definition: basics.h:448
void check_nullspace_ok()
Definition: sim_utils.cc:202
SF::vector< stimulus > stimuli
the electrical stimuli
Definition: electrics.h:254
const char * name
timer name
Definition: timer_utils.h:53
void update_cwd()
save the current working directory to curdir so that we can switch back to it if needed.
Definition: sim_utils.cc:478
bool phys_defined(int physreg)
function to check if certain physics are defined
File descriptor struct.
Definition: basics.h:133
void redistribute_elements(meshdata< T, S > &mesh, meshdata< T, S > &sendbuff, vector< T > &part)
Redistribute the element data of a parallel mesh among the ranks based on a partitioning.
int set_dir(IO_t dest)
Definition: sim_utils.cc:483
virtual void output_timings()
Definition: physics_types.h:76
std::string get_basename(const std::string &path)
Definition: basics.cc:61
std::string name
the mesh name
Definition: SF_container.h:395
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
Top-level header of physics module.
void log_msg(FILE_SPEC out, int level, unsigned char flag, const char *fmt,...)
Definition: basics.cc:72
#define NONL
Definition: basics.h:312
void ignore_extracellular_stim(Stimulus *st, int ns, int ignore)
Definition: sim_utils.cc:154
#define MONODOMAIN
Definition: sim_utils.h:143
bool mesh_is_registered(const mesh_t gt)
check wheter a SF mesh is set
Definition: sf_interface.cc:59
stim_t type
type of stimulus
Definition: stimulate.h:136
bool have_permutation(const int mesh_id, const int perm_id, const int dpn)
short get_mesh_dim(mesh_t id)
get (lowest) dimension of the mesh used in the experiment
Definition: sim_utils.cc:1299
SF::vector< async_io_item > async_IOs
Definition: sim_utils.h:306
A vector storing arbitrary data.
Definition: SF_vector.h:42
size_t l_numelem
local number of elements
Definition: SF_container.h:387
#define STM_IGNORE_PSEUDO_BIDM
Definition: stimulate.h:61
#define Extracellular_V_OL
Definition: stimulate.h:42
std::map< datavec_t, sf_vec * > datavec_reg
important solution vectors from different physics
Definition: main.cc:54
virtual void compute_step()=0
#define Transmembrane_I
Definition: stimulate.h:37
void f_close(FILE_SPEC &f)
Close a FILE_SPEC.
Definition: basics.cc:162
std::vector< base_timer * > timers
vector containing individual timers
Definition: timer_utils.h:84
limpet::MULTI_IF * miif
Definition: ionics.h:66
std::string timer_unit(const int timer_id)
figure out units of a signal linked to a given timer
Definition: electrics.cc:733
void initialize_neq_timer(const std::vector< double > &itrig, double idur, int ID, const char *iname, const char *poolname=nullptr)
Definition: timer_utils.cc:63
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:1447
std::map< int, std::string > units
Definition: stimulate.cc:41
#define Extracellular_V
Definition: stimulate.h:39
int load_ionic_module(const char *)
const vector< T > & algebraic_nodes() const
Getter function for the local indices forming the local algebraic node set.
vector< S > she
sheet direction
Definition: SF_container.h:409
datavec_t
Enum used to adress the different data vectors stored in the data registry.
Definition: physics_types.h:99
SF::vector< sync_io_item > sync_IOs
Definition: sim_utils.h:305
double strength
strength of stimulus
Definition: stimulate.h:92
stim_protocol ptcl
applied stimulation protocol used
Definition: stimulate.h:166
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
Functor class applying a submesh renumbering.
Definition: SF_numbering.h:67
Basic utility structs and functions, mostly IO related.
stim_physics phys
physics of stimulus
Definition: stimulate.h:167
#define Extracellular_Ground
Definition: stimulate.h:40
size_t root_write(FILE *fd, const vector< V > &vec, MPI_Comm comm)
Write vector data binary to disk.
void unites_y(const char *a)
Definition: IGBheader.h:348
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
sf_vec * get_data(datavec_t d)
Retrieve a petsc data vector from the data registry.
Definition: sim_utils.cc:880
#define IGB_VEC4_f
Definition: IGBheader.h:73
int get_size(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:290
int IO_id
pointer to data registered for output
Definition: sim_utils.h:282
Simulator-level utility execution control functions.
double start
initial time (nonzero when restarting)
Definition: timer_utils.h:79
SF::meshdata< mesh_int_t, mesh_real_t > sf_mesh
Definition: sf_interface.h:47
centralize time managment and output triggering
Definition: timer_utils.h:73
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:276
SF::scatter_registry scatter_reg
Registry for the different scatter objects.
Definition: main.cc:44
int postproc_recover_phie()
Definition: electrics.cc:1936
T local_nodal_to_local_petsc(const meshdata< T, S > &mesh, int rank, T local_nodal)
virtual T gsize() const =0
virtual void initialize()=0
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:135
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:844
Container for a PETSc VecScatter.