openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
MULTI_ION_IF.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 
54 // TODO remove the "get_sv_address()" and "get_sv_size()" calls by moving
55 // the functions that use them to the IonIf class
56 
57 #include "MULTI_ION_IF.h"
58 #include "limpet_types.h"
59 #include <stdio.h>
60 #ifdef HAS_CUDA_MODEL
61 #include <cuda_runtime.h>
62 #endif
63 #if defined HAS_ROCM_MODEL && defined __HIP__
64 #include <hip/hip_runtime.h>
65 #endif
66 
67 #include "SF_init.h" // for SF::init_xxx()
68 #include "petsc_utils.h" // TODO for EXIT
69 
70 namespace limpet {
71 
72 using ::opencarp::base_timer;
83 using ::opencarp::Salt_list;
88 #ifdef USE_FMEM_WRAPPER
90 #endif
91 
92 /*
93  PROTOTYPES
94  */
95 void initialize_params_MIIF(MULTI_IF *pMIIF);
96 void initialize_ionic_IF(MULTI_IF *pMIIF);
97 void CreateIIFLocalNodeLsts(MULTI_IF *pMIIF);
98 void CreateIIFNodeLsts_(int, char *, int **, int ***, int);
99 void alloc_MIIF(MULTI_IF *pMIIF);
100 void initializeIMPData(MULTI_IF *pMIIF);
101 void freeIMPData(MULTI_IF *pMIIF);
102 void allocate_shared_data(MULTI_IF *);
103 int getGlobalNodalIndex(IonIfBase& pIF, int relIdx);
105 bool isIMPdata(const char *);
106 
107 static MULTI_IF *gMIIF_Error_Recovery;
108 static int current_IIF_index;
109 static float current_time = 0;
110 static float start_time = 0;
111 
112 static void prepare_error_recovery(MULTI_IF *MIIF, int IIF_index, float time) {
113  gMIIF_Error_Recovery = MIIF;
114  current_IIF_index = IIF_index;
115  current_time = start_time+time;
116 }
117 
118 int current_global_node(int local_node) {
119  return gMIIF_Error_Recovery->NodeLists[current_IIF_index][local_node];
120 }
121 
123  return current_time;
124 }
125 
126 static int g_print_bounds_exceeded = 1;
127 
129  return g_print_bounds_exceeded;
130 }
131 
133  int oldval = g_print_bounds_exceeded;
134 
135  g_print_bounds_exceeded = newval;
136  return oldval;
137 }
138 
142  CreateIIFNodeLsts_(pMIIF->N_IIF, pMIIF->IIFmask, &pMIIF->N_Nodes,
143  &pMIIF->NodeLists, pMIIF->numNode);
144 }
145 
157 void CreateIIFNodeLsts_(int N_IIF, IIF_Mask_t *IIF_Mask, int **N_Nodes,
158  int ***NodeLists, int numNode)
159 {
160  /* create node lists for each IIF
161  NodeNum: store number of nodes for each IIF
162  NodeLst: store node indices of each IIF
163  */
164  int *NodeNum = static_cast<int *>(calloc(N_IIF, sizeof(int)));
165  int **NodeLst = static_cast<int **>(calloc(N_IIF, sizeof(int *)));
166 
167  // determine number of nodes of each IIF
168  for (int i = 0; i < numNode; i++)
169  NodeNum[static_cast<int>(IIF_Mask[i])]++;
170 
171  /* now we know the number of nodes per IIF and store the indices of
172  nodes belonging to a particular IIF into a list.
173  */
174  for (int i = 0; i < N_IIF; i++) {
175  if (NodeNum[i] > 0)
176  NodeLst[i] = static_cast<int *>(malloc(sizeof(int *) * NodeNum[i]));
177  else
178  NodeLst[i] = NULL;
179  }
180 
181  // just store indices relative to low
182  for (int j = 0; j < N_IIF; j++) {
183  int hcount = 0;
184  for (int i = 0; i < numNode; i++)
185  if (IIF_Mask[i] == j)
186  NodeLst[j][hcount++] = i;
187  }
188 
189  for(int j=0; j<N_IIF; j++)
190  std::sort(NodeLst[j], NodeLst[j]+NodeNum[j]);
191 
192  *N_Nodes = NodeNum;
193  *NodeLists = NodeLst;
194 
195 } // CreateIIFNodeLsts_
196 
197 
205 int getGlobalNodalIndex(IonIfBase& pIF, int rIdx) {
206  // get pointer to parent MIIF
207  MULTI_IF *pMIIF = (MULTI_IF *) pIF.parent(); // TODO check?!?!?!?!?!
208 
209  if (pIF.get_type().is_plugin())
210  pMIIF = (MULTI_IF *) pIF.parent()->parent();
211 
212  return pMIIF->NodeLists[pIF.miifIdx][rIdx];
213 }
214 
221 void alloc_MIIF(MULTI_IF *pMIIF) {
222  int N_IIF = pMIIF->N_IIF;
223 
224  pMIIF->contiguous = static_cast<bool *>(calloc(N_IIF, sizeof(bool)));
225 #if 0
226  pMIIF->ldata = (GlobalData_t***) build_matrix_ns<GlobalData_t>(N_IIF, NUM_IMP_DATA_TYPES, sizeof(GlobalData_t **), Target::CPU);
227 #else
228  pMIIF->ldata = allocate_on_target<GlobalData_t**>(Target::CPU, N_IIF);
229  for (std::size_t i = 0; i < pMIIF->N_IIF; i++) {
230  pMIIF->ldata[i] = allocate_on_target<GlobalData_t*>(pMIIF->iontypes[i].get().select_target(pMIIF->targets[i]), NUM_IMP_DATA_TYPES);
231  }
232 #endif
233 
234  for (int i = 0; i < pMIIF->N_IIF; i++) {
235  // Create an IonIf object on the chosen target
236  pMIIF->IIF.push_back(pMIIF->iontypes[i].get().make_ion_if(pMIIF->targets[i], pMIIF->N_Nodes[i], pMIIF->plugtypes[i]));
237  // this is questionable. A IIF does not have a parent, parent should be NULL
238  // It might make sense in some situtation to be able to refer to the MIIF
239  // structure, but pIF->parent is of type IonIf and not MULTI_IF!!!!
240  pMIIF->IIF[i]->set_parent((IonIfBase *)pMIIF);
241 
242  // set unique miifIdx for imp and plugins to enable refering back
243  // to global entities such as global node index from within imp
244  pMIIF->IIF[i]->for_each([&](IonIfBase& IF) { IF.miifIdx = i; });
245  }
246 } // alloc_MIIF
247 
253  for (auto& imp : pMIIF->IIF) {
254  imp->initialize_params();
255  }
256 }
257 
262  pMIIF->getRealData(); // needed for initializing Vm, etc.
263  for (int i = 0; i < pMIIF->IIF.size(); i++)
264  pMIIF->IIF[i]->initialize(pMIIF->dt, pMIIF->ldata[i]);
265 
266  pMIIF->releaseRealDataDuringInit();
267 }
268 
275  // we need to copy all global data which has been set, even if it is not
276  // listed as moddat. This is the case with Vm which may not be modified by
277  // LIMPET but is initialized by LIMPET
278  unsigned int *moddat = static_cast<unsigned int *>(calloc(this->N_IIF, sizeof(unsigned int)));
279 
280  for (int i = 0; i < this->N_IIF; i++) {
281  moddat[i] = this->IIF[i]->get_moddat();
282  this->IIF[i]->set_moddat(moddat[i] | this->IIF[i]->get_reqdat()); // add in the reqdat
283  }
284  this->releaseRealData();
285  for (int i = 0; i < this->N_IIF; i++) {
286  this->IIF[i]->set_moddat(moddat[i]); // remove the reqdat
287  }
288  free(moddat);
289 }
290 
296 {
298  alloc_MIIF(this);
300 
301  memset(&svd, 0, sizeof(SV_DUMP) );
302  svd.active = 0;
303  svd.intv = 1.0; // 1. ms by default
304 }
305 
306 #define FILENAME_BUF 1024
307 
309 int int_cmp(const void *a, const void *b) {
310  return *(int *)a - *(int *)b;
311 }
312 
340 void open_trace(MULTI_IF *MIIF, int n_traceNodes, int *traceNodes, int *label, opencarp::sf_mesh* imesh)
341 {
342  if (!n_traceNodes) return;
343 
344  if (n_traceNodes > 1000)
345  log_msg(0, 4, 0, "%s warning: %d trace nodes may impact performance", __func__, n_traceNodes);
346 
347  MIIF->trace_info = (Trace_Info *)calloc(n_traceNodes+1, sizeof(Trace_Info));
348  MIIF->trace_info[n_traceNodes].region = -1; // end of list marker
349 
350  // bijective index mapping between set A (local petsc indexing) and set B (local nodal indexing)
352  if(imesh)
353  SF::local_petsc_to_nodal_mapping(*imesh, petsc2nod);
354 
355  for (int iTrace = 0; iTrace < n_traceNodes; iTrace++) {
356  mesh_int_t lnode = imesh ? imesh->pl.localize(traceNodes[iTrace]) : traceNodes[iTrace];
357 
358  // here we check if lnode is in set B (i.e. local nodal indexing)
359  if(imesh) {
360  // here we check if lnode is in set B (i.e. local nodal indexing)
361  if(petsc2nod.in_b(lnode) == false) continue;
362 
363  // lnode is local nodal, we map it to local petsc
364  lnode = petsc2nod.backward_map(lnode);
365  }
366 
367  for (int iRegion = 0; iRegion < MIIF->N_IIF; iRegion++)
368  {
369  auto target = static_cast<int *>(bsearch(&lnode, MIIF->NodeLists[iRegion],
370  MIIF->N_Nodes[iRegion], sizeof(int), int_cmp));
371 
372  if (target != NULL ) {
373  int idx = target - MIIF->NodeLists[iRegion];
374  Trace_Info *trace_info = MIIF->trace_info+iTrace;
375  trace_info->found = true;
376  trace_info->region = iRegion;
377  trace_info->node_idx = idx;
378  }
379  }
380  }
381 
382  for (int iTrace = 0; iTrace < n_traceNodes; iTrace++) {
383  if(get_global(int(MIIF->trace_info[iTrace].found), MPI_SUM) == 0) {
384  MIIF->trace_info[iTrace].ignored = true;
385  log_msg(0,4,0, "trace node %d not found", traceNodes[iTrace]);
386  continue;
387  }
388 
389  char traceName[FILENAME_BUF];
390  snprintf(traceName, sizeof traceName, "Trace_%d.dat", label ? label[iTrace] : traceNodes[iTrace]);
391  MIIF->trace_info[iTrace].file = f_open(traceName, "w");
392  }
393 } // open_trace
394 
409 void dump_trace(MULTI_IF *MIIF, limpet::Real time) {
410  if (MIIF->trace_info == NULL)
411  return;
412 
413 // struct exception_type e;
414 
415 #define MAX_TRACE_LINE_LEN 8196
416 
417  char trace_buf[MAX_TRACE_LINE_LEN] = {0};
418 
419  std::vector<IonIfBase*>& IIF = MIIF->IIF;
420  FILE *fs;
421 
422  for (int iTrace = 0; MIIF->trace_info[iTrace].region >= 0; iTrace++) {
423  Trace_Info *ctrace = MIIF->trace_info+iTrace;
424 
425  if (ctrace->ignored)
426  continue;
427 
428  if (ctrace->found) {
429  fs = fmemopen_(trace_buf, MAX_TRACE_LINE_LEN, "w");
430  fprintf(fs, "%4.10f\t", time);
431  if (IIF[ctrace->region]->get_type().has_trace()) {
432  IIF[ctrace->region]->get_type().trace(*IIF[ctrace->region], ctrace->node_idx,
433  fs, MIIF->ldata[ctrace->region]);
434  }
435  for (auto& plugin : IIF[ctrace->region]->plugins()) {
436  if (plugin->get_type().has_trace()) {
437  fprintf(fs, "\t");
438  plugin->get_type().trace(
439  *plugin, ctrace->node_idx,
440  fs, MIIF->ldata[ctrace->region]);
441  }
442  }
443  fprintf(fs, "\n");
444  fclose(fs);
445 
446  fprintf(ctrace->file->fd, "%s", trace_buf);
447  }
448  fflush(ctrace->file->fd);
449  }
450 } // dump_trace
451 
452 void close_trace(MULTI_IF *MIIF) {
453  if (!MIIF->trace_info) return;
454 
455  for (int iTrace = 0; MIIF->trace_info[iTrace].region >= 0; iTrace++)
456  f_close(MIIF->trace_info[iTrace].file);
457 
458  free(MIIF->trace_info);
459  MIIF->trace_info = NULL;
460 }
461 
468 void MULTI_IF::initialize_currents(double idt, int subDt) {
469  numSubDt = subDt;
470  dt = idt/this->numSubDt;
471 
472  allocate_shared_data(this);
473  initializeIMPData(this);
474  initialize_ionic_IF(this);
475 }
476 
483 void MULTI_IF::dump_luts_MIIF(bool zipped) {
484  // if one of the IFs in MIIF does not live in partion 0
485  // we won't get a LUT dumped.
486  if (get_rank()) return;
487 
488  std::vector<IonIfBase*>& pIF = this->IIF;
489  for (auto& IF : pIF) {
490  int ndmps = IF->dump_luts(zipped);
491  if (ndmps < IF->tables().size()) {
492  log_msg(logger, 4, 0, "LUT dump error %s: only %d out of %d LUTs dumped.\n",
493  IF->get_type().get_name().c_str(), ndmps, IF->tables().size());
494  }
495  for (auto& plugin : IF->plugins()) {
496  ndmps = plugin->dump_luts(zipped);
497  if (ndmps < plugin->tables().size()) {
498  log_msg(logger, 4, 0, "LUT dump error %s: only %d out of %d LUTs dumped.\n",
499  plugin->get_type().get_name().c_str(), ndmps, IF->tables().size());
500  }
501  }
502  }
503 }
504 
510  if (get_rank()) return;
511 
512  // close files
513  for (int i = 0; i < svd.n; i++)
514 #ifndef USE_HDF5
515  if (svd.hdls[i] != NULL)
516 #endif // ifndef USE_HDF5
517  f_close(svd.hdls[i]);
518 }
519 
530 char *get_sv(void *tab, int offset, int n, int svSize, int size) {
531  char *buf = static_cast<char *>(malloc(n*size));
532  char *bp = buf;
533  char *p = static_cast<char *>(tab) + offset;
534 
535  for (int i = 0; i < n; i++) {
536  memcpy(bp, p, size);
537  bp += size;
538  p += svSize;
539  }
540  return buf;
541 }
542 
550 int MULTI_IF::dump_svs(base_timer *iot) {
551  int nwr = 0;
552 
553  if (iot->triggered) {
554  for (int i = 0; i < svd.n; i++) {
555  nwr = 0;
556 
557  FILE* fd = svd.hdls[i] ? svd.hdls[i]->fd : NULL;
558  char *buf = get_sv(svd.svtab[i], svd.offset[i], svd.num[i], svd.svsize[i], svd.size[i]);
559  nwr += SF::root_write<char>(fd, (char*) buf, svd.size[i]*svd.num[i], PETSC_COMM_WORLD);
560  free(buf);
561  }
562  svd.nwr += nwr;
563  svd.n_dumps++;
564  }
565 
566  return nwr;
567 }
568 
574 #ifdef HAS_GPU_MODEL
575 __global__
576 void update_vm(int start, int end, double *vm, double *ion, double dt)
577 {
578  int i = blockIdx.x*blockDim.x + threadIdx.x;
579  if (i < end)
580  vm[i] = vm[i] + (ion[i] * (-dt));
581 }
582 #endif
583 
584 // * compute what has to be computed, current or otherwise
585 void MULTI_IF::compute_ionic_current(bool flag_send, bool flag_receive)
586 {
587  gdata[Iion]->set(0.0);
588  if (flag_send == 1) {
589  this->getRealData();
590  }
591 
592  for (int j = 0; j < this->numSubDt; j++)
593  {
594  for (int i = 0; i < this->N_IIF; i++)
595  {
596  IonIfBase* pIF = this->IIF[i];
597  if (!this->N_Nodes[i]) continue;
598 
599  prepare_error_recovery(this, i, pIF->get_tstp().cnt * this->dt);
600  update_ts(&pIF->get_tstp());
601 
602  int current = 0;
603  do {
604  try {
605  pIF->compute(current, this->N_Nodes[i], this->ldata[i]);
606  current = this->N_Nodes[i];
607  }
608  catch(int e) { //FIXME
609  if(e==-1) {
610  // CHECK this->NodeLists[i][e.node], e.node);
611  fprintf(stderr, "LIMPET compute fail in %s at node %d (local %d)! Aborting!\n",
612  this->name.c_str(), this->NodeLists[i][j], j);
613  exit(1);
614  }
615  current++;
616  }
617  } while (current < this->N_Nodes[i]);
618 
619  for (auto& plugin : pIF->plugins()) {
620  current = 0;
621  update_ts(&plugin->get_tstp());
622 
623  do {
624  try {
625  plugin->compute(current, this->N_Nodes[i], this->ldata[i]);
626  current = this->N_Nodes[i];
627 
628  if (plugin == nullptr)
629  throw -1;
630  }
631  catch(int e) { //FIXME
632  if(e==-1) {
633  // CHECK this->NodeLists[i][e.node], e.node);
634  fprintf(stderr, "LIMPET plugin fail in %s at node %d (local %d)! Aborting!\n",
635  this->name.c_str(), this->NodeLists[i][j], j);
636  exit(1);
637  }
638  current++;
639  }
640  } while (current < this->N_Nodes[i]);
641  }
642  }
643 
644  #ifdef HAS_GPU_MODEL
645  //TODO: extUpdateVm can be true outside of the bench executable ! But this should
646  // only run in bench! This should be done for each model according to its target
647  if(!extUpdateVm && is_gpu(this->targets[0])) {
648 
649  #if defined __CUDA__ || defined __HIP__
650  update_vm<<<(this->N_Nodes[0]/64)+1,64>>>(0, this->N_Nodes[0], this->ldata[0][Vm], this->ldata[0][Iion], this->dt);
651  #ifdef __CUDA__
652  cudaDeviceSynchronize();
653  #elif defined __HIP__
654  hipDeviceSynchronize();
655 #endif
656  #else
657  fprintf(stderr, "GPU/CUDA not found");
658  #endif
659  }
660  #endif
661 
662  if (flag_receive == 1) {
663  this->releaseRealData();
664  }
665 
666  // TODO: Target should be checked for each region
667  if(!extUpdateVm && !is_gpu(this->targets[0]))
668  gdata[Vm]->add_scaled(*gdata[Iion], SF_real(-dt));
669  }
670 }
671 
673  assert(!this->doppel);
674 
675  // Free data first since this function needs to access memory in the IIFs
676  freeIMPData(this);
677  for (auto& pIF : this->IIF) {
678  pIF->get_type().destroy_ion_if(pIF);
679  }
680 }
681 
689  for (int i = 0; i < pMIIF->N_IIF; i++) {
690  if (!pMIIF->N_Nodes[i])
691  continue;
692 
693  // if memory contiguous, we don't need to copy data, merely pass pointers
694  pMIIF->contiguous[i] = pMIIF->N_Nodes[i]-1 ==
695  pMIIF->NodeLists[i][pMIIF->N_Nodes[i]-1]-pMIIF->NodeLists[i][0];
696 
697  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
698  // check if it is used by the IMPs
699  if (!USED_DAT(pMIIF->IIF[i], imp_data_flag[j]) )
700  continue;
701 
702  // allocate mem buffer only if not contiguous or if data is on GPU
703  if (!pMIIF->contiguous[i] || is_gpu(pMIIF->IIF[i]->get_target()))
704  pMIIF->ldata[i][j] =
705  allocate_on_target<GlobalData_t>(pMIIF->IIF[i]->get_target(),
706  pMIIF->N_Nodes[i]);
707  // check if data supplied
708  if (pMIIF->gdata[j] == NULL) {
709  log_msg(pMIIF->logger, 5, LOCAL, "IMP data type %s not supplied for region %d", imp_data_names[j], i);
710  exit(1);
711  }
712  }
713  }
714 } // initializeIMPData
715 
725  // set rdata in the parent vector to point to the
726  // local data to be computed on the local processor
727  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
728  SF_real* rdata;
729  if (this->gdata[j] != NULL) {
730  rdata = this->gdata[j]->ptr();
731  this->procdata[j] = rdata;
732  }
733  else continue;
734 
735  // now map local data so that each ionic model can use it
736  // if noncontiguous in the parent vector, copy to contiguous local vector
737  for (int i = 0; i < this->N_IIF; i++) {
738  if (!this->N_Nodes[i] || !USED_DAT(this->IIF[i], imp_data_flag[j]) )
739  continue;
740 
741  if (this->contiguous[i] && !is_gpu(this->IIF[i]->get_target()))
742  this->ldata[i][j] = static_cast<GlobalData_t *>(rdata) + this->NodeLists[i][0];
743  else {
744  const int* ip = this->NodeLists[i];
745  for (int k = 0; k < this->N_Nodes[i]; k++)
746  this->ldata[i][j][k] = rdata[ip[k]];
747  }
748  }
749  }
750 } // getRealData
751 
761 
762  for (int i = 0; i < this->N_IIF; i++) {
763  if (!this->N_Nodes[i]) continue;
764 
765  // Every external variable used by the model must be copied on GPU
766  if (is_gpu(this->IIF[i]->get_target())) {
767  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
768  if (!USED_DAT(this->IIF[i], imp_data_flag[j]) )
769  continue;
770  if (this->contiguous[i])
771  memcpy(&this->procdata[j][this->NodeLists[i][0]], this->ldata[i][j], sizeof(this->ldata[i][j][0])*this->N_Nodes[i]);
772  else {
773  for (int k = 0; k < this->N_Nodes[i]; k++)
774  this->procdata[j][this->NodeLists[i][k]] = this->ldata[i][j][k];
775  }
776  }
777  }
778  else {
779  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
780  if (this->IIF[i]->get_moddat() & imp_data_flag[j]) {
781  if (!this->contiguous[i]) {
782  for (int k = 0; k < this->N_Nodes[i]; k++)
783  this->procdata[j][this->NodeLists[i][k]] = this->ldata[i][j][k];
784  }
785  }
786  }
787  }
788  }
789 
790  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++)
791  if (this->gdata[j] != NULL)
792  this->gdata[j]->release_ptr(this->procdata[j]);
793 }
794 
795 // * free the local storage vectors
796 void freeIMPData(MULTI_IF *pMIIF) {
797  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
798  if (pMIIF->gdata[j] != NULL) {
799  for (int i = 0; i < pMIIF->N_IIF; i++)
800  if ((!pMIIF->contiguous[i] || is_gpu(pMIIF->IIF[i]->get_target())) && USED_DAT(pMIIF->IIF[i], imp_data_flag[j]) && pMIIF->N_Nodes[i] > 0) {
801  deallocate_on_target<GlobalData_t>(pMIIF->IIF[i]->get_target(),
802  pMIIF->ldata[i][j]);
803  }
804  delete pMIIF->gdata[j];
805  }
806  }
807  free(pMIIF->contiguous);
809 }
810 
820 int get_plug_flag(char *plgstr, int *out_num_plugins, IonTypeList& out_plugins) {
821  if (plgstr == NULL) {
822  *out_num_plugins = 0;
823  out_plugins.clear();
824  return 1;
825  }
826 
827  // maximum number of plugins is 200 apparently
828  char *plugspec = dupstr(plgstr);
829  char *saveptr;
830  char *token = tokstr_r(plugspec, ":", &saveptr);
831 
832  *out_num_plugins = 0;
833  while (token) {
834  IonType* type = get_ion_type(std::string(token));
835  if (type == NULL) {
836  free(plugspec);
837  return 0;
838  } else {
839  out_plugins.push_back(*type);
840  }
841  token = tokstr_r(NULL, ":", &saveptr);
842  }
843 
844  *out_num_plugins = out_plugins.size();
845 
846  free(plugspec);
847  return 1;
848 } // get_plug_flag
849 
868 int MULTI_IF::adjust_MIIF_variables(const char* variable,
869  const SF::vector<int> & indices,
870  const SF::vector<double> & values)
871 {
872  int num_changed = 0;
873 
874  // determine if we're dealing with an external variable or a state variable.
875  if (index(variable, '.') == NULL) {
876  // external variable, the easy case.
877  // which external variable?
878  int data_id = -1;
879  for (int ii = 0; ii < NUM_IMP_DATA_TYPES; ii++) {
880  if (strcmp(variable, imp_data_names[ii]) == 0) {
881  data_id = ii;
882  break;
883  }
884  }
885  if (data_id == -1) {
886  log_msg(logger, 5, FLUSH, "Error! No external variable named %s in this build of openCARP.", variable);
887  exit(EXIT_FAILURE);
888  }
889 
890  // make sure that external variable is being used.
891  if (this->gdata[data_id] == NULL) {
892  log_msg(logger, 4, 0,
893  "External variable %s is not being used this processor.\n"
894  "Is this really what you meant to do?\n"
895  "Perhaps you don't have the correct Ionic models selected.",
896  imp_data_names[data_id]);
897  }
898  else {
899  // fill in the external variable with the data from the file.
900  SF_real *raw_data = this->gdata[data_id]->ptr();
901 
902  for (size_t i = 0; i < indices.size(); i++) {
903  raw_data[indices[i]] = values[i];
904  num_changed++;
905  }
906 
907  this->gdata[data_id]->release_ptr(raw_data);
908  }
909  }
910  else {
911  // state variable.
912  // extract the IMP name and state variable name from the string.
913  char *saveptr;
914  char *my_variable = dupstr(variable);
915  char *IIF_name = tokstr_r(my_variable, ".", &saveptr);
916  char *sv_name = tokstr_r(NULL, ".", &saveptr);
917 
918  IonType* type = get_ion_type(std::string(IIF_name));
919  if (type == NULL) {
920  log_msg(logger, 5, 0, "%s error: %s is not a valid IMP name.", __func__, IIF_name);
921  exit(EXIT_FAILURE);
922  }
923 
924  int sv_offset;
925  int sv_size;
926  SVgetfcn sv_get = type->get_sv_offset(sv_name, &sv_offset, &sv_size);
927  if (sv_get == NULL) {
928  log_msg(logger, 5, 0, "%s error: %s is not a valid state variable for the %s model.",
929  __func__, sv_name, IIF_name);
930  exit(EXIT_FAILURE);
931  }
932  SVputfcn sv_put = getPutSV(sv_get);
933 
934  // Ok, go through the ionic models
935  for (int i_iif = 0; i_iif < this->N_IIF; i_iif++) {
936  IonIfBase *lIIF = NULL;
937  if (this->IIF[i_iif]->get_type() == *type) {
938  lIIF = this->IIF[i_iif];
939  }
940  else {
941  // Go through the plugins
942  for (auto& plugin : this->IIF[i_iif]->plugins()) {
943  if (plugin->get_type() == *type) {
944  lIIF = plugin;
945  break;
946  }
947  }
948  }
949 
950  if (lIIF == NULL) {
951  continue;
952  }
953 
954  for (size_t ii = 0; ii < indices.size(); ii++) {
955  GlobalData_t file_value = values[ii];
956 
957  int *target = static_cast<int*>(bsearch(&indices[ii], this->NodeLists[i_iif],
958  this->N_Nodes[i_iif], sizeof(int), int_cmp));
959  if (target) {
960  // We found a point to adjust! Do the adjustment.
961  num_changed++;
962  sv_put(*lIIF, (int)(target - this->NodeLists[i_iif]), sv_offset, file_value);
963  }
964  }
965 
966  }
967  free(my_variable);
968  }
969 
970  MPI_Allreduce(MPI_IN_PLACE, &num_changed, 1, MPI_INT, MPI_SUM, PETSC_COMM_WORLD);
971 
972  return num_changed;
973 } // adjust_MIIF_variables
974 
992 int determine_write_ranges(int N, size_t *offset, size_t bufsize, int **ranges) {
993  int nitems;
994 
995  if (!get_rank() ) {
996  Salt_list r; STRUCT_ZERO(r); r.chunk = 10;
997  int temp = 0;
998  SLIST_APPEND(&r, temp);
999  long loff = offset[0];
1000  for (int i = 0; i < N; i++)
1001  if (offset[i]-loff > bufsize) {
1002  SLIST_APPEND(&r, i);
1003  loff = offset[i];
1004  }
1005 
1006  nitems = r.nitems;
1007  SLIST_APPEND(&r, N);
1008  *ranges = (int *)r.data;
1009  }
1010 
1011  MPI_Bcast(&nitems, 1, MPI_INT, 0, PETSC_COMM_WORLD);
1012  if (get_rank() )
1013  *ranges = static_cast<int *>(malloc( (nitems+1)*sizeof(int) ));
1014  MPI_Bcast(*ranges, nitems+1, MPI_INT, 0, PETSC_COMM_WORLD);
1015  return nitems;
1016 } // determine_write_ranges
1017 
1056 void MULTI_IF::dump_state(char *fname, float simtime, mesh_t gid,
1057  bool append, unsigned int revision)
1058 {
1059  float t0, t1;
1060  t0 = get_time();
1061 
1062  FILE_SPEC out = NULL;
1063  int rank = get_rank();
1064  int miif_node_gsize = get_global(this->numNode, MPI_SUM);
1065  int error = 0;
1066 
1067  log_msg(logger, 0, 0, "Saving state at time %f in file: %s", simtime, fname);
1068 
1069  if (rank == 0) {
1070  out = f_open(fname, append ? "a" : "w");
1071 
1072  if(out) {
1073  fseek(out->fd, 0, SEEK_END); // NOP call to sync disk data for switch to write mode
1074  write_bin_string(out, Magic_MIIF_ID);
1075 
1076  fwrite(&MIIF_Format, sizeof(unsigned int), 1, out->fd);
1077  fwrite(&revision, sizeof(unsigned int), 1, out->fd);
1078  time_t tm = time(NULL);
1079  fwrite(&tm, sizeof(time_t), 1, out->fd);
1080  fwrite(&simtime, sizeof(float), 1, out->fd);
1081  fwrite(&miif_node_gsize, sizeof(int), 1, out->fd);
1082 
1083  int num_gdata = 0;
1084  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++)
1085  if (this->gdata[i]) num_gdata++;
1086 
1087  fwrite(&num_gdata, sizeof(int), 1, out->fd);
1088  }
1089  else
1090  error++;
1091  }
1092 
1093 // #define CHATTY
1094 
1095  if(get_global(error, MPI_SUM))
1096  EXIT(EXIT_FAILURE);
1097 
1098  FILE* fd = rank == 0 ? out->fd : NULL;
1099  sf_vec* outVec;
1100 
1101  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++) {
1102  if (this->gdata[i]) {
1103  if(rank == 0) {
1104  write_bin_string(out, imp_data_names[i]);
1105 #ifdef CHATTY
1106  log_msg(NULL, 0, 0, "\tDumping %s at %d", imp_data_names[i], ftell(fd) );
1107 #endif
1108  }
1109 
1110  // If we have a proper set up parallel context, we pass a gid != unset_msh. Then
1111  // we can query petsc-to-canonical scattering. If we are in a simpler context
1112  // like bench, we assume that no scattering is needed and just shallow-copy the vector.
1113  if(gid != unset_msh) {
1114  SF::init_vector(&outVec, this->gdata[Vm]);
1116  assert(sc != NULL);
1117 
1118  sc->forward(*this->gdata[i], *outVec);
1119  outVec->write_binary<SF_real>(fd);
1120  delete outVec;
1121  }
1122  else {
1123  this->gdata[i]->write_binary<SF_real>(fd);
1124  }
1125  }
1126  }
1127 
1128  // output IMP region info
1129  // #regions and the IM and plugins for each region
1130  // determine memory requirements for each region
1131  int* imp_mem = new int[this->N_IIF];
1132  size_t *offset = NULL;
1133  long filepos = 0.;
1134 
1135  if (rank == 0) {
1136 #ifdef CHATTY
1137  log_msg(NULL, 0, 0, "\tDumping IMP sizes at %d", ftell(fd) );
1138 #endif
1139 
1140  int max = 1;
1141  fwrite(&this->N_IIF, sizeof(int), 1, fd); // #IIF's
1142 
1143  for (int i = 0; i < this->N_IIF; i++) {
1144  IonIfBase *imp = this->IIF[i];
1145  write_bin_string(out, imp->get_type().get_name().c_str());
1146 
1147  std::size_t sv_size = imp->get_sv_size();
1148  fwrite(&sv_size, sizeof(int), 1, fd);
1149  unsigned long n_plugins = (unsigned long) imp->plugins().size();
1150  fwrite(&n_plugins, sizeof(int), 1, fd);
1151  imp_mem[i] = imp->get_sv_size();
1152 
1153  for (auto& plug : imp->plugins()) {
1154  write_bin_string(out, plug->get_type().get_name().c_str());
1155 
1156  fwrite(&sv_size, sizeof(int), 1, fd);
1157  imp_mem[i] += plug->get_sv_size();
1158  }
1159  }
1160  filepos = ftell(fd);
1161  }
1162 
1163  MPI_Bcast(imp_mem, this->N_IIF, MPI_INT, 0, PETSC_COMM_WORLD);
1164 
1165  // If we have a proper set up parallel context, we pass a gid != unset_msh. Then
1166  // we can query a mesh and its parallel numbering. If we are in a simpler context
1167  // like bench, we build the numbering ourselves. Then, we assume that no
1168  // renumbering has taken place.
1169  SF::vector<int> loc2canon;
1170  if(gid != unset_msh) {
1171  const sf_mesh & mesh = get_mesh(gid);
1172  const SF::vector<mesh_int_t> & canon_nbr = mesh.get_numbering(SF::NBR_SUBMESH);
1173  const SF::vector<mesh_int_t> & alg_idx = mesh.pl.algebraic_nodes();
1174  loc2canon.resize(alg_idx.size());
1175 
1176  for(size_t i=0; i<alg_idx.size(); i++) loc2canon[i] = canon_nbr[alg_idx[i]];
1177  }
1178  else {
1179  int rank = get_rank();
1180  SF::vector<int> layout;
1181  layout_from_count(this->numNode, layout, PETSC_COMM_WORLD);
1182  loc2canon.resize(this->numNode);
1183 
1184  for (int i = 0; i < this->numNode; i++) loc2canon[i] = layout[rank] + i;
1185  }
1186 
1187  // write out IIF_Mask canonically ordered
1188  int *canord = new int[this->numNode];
1189  for (int i = 0; i < this->numNode; i++)
1190  canord[i] = loc2canon[i];
1191 
1192 #ifdef CHATTY
1193  log_msg(NULL, 0, 0, "\tDumping IMP masks at %d", filepos);
1194 #endif // ifdef CHATTY
1195 
1196  SF::root_write_ordered<int, IIF_Mask_t>(fd, canord, IIFmask, numNode, PETSC_COMM_WORLD);
1197 
1198  // for every local node, we compute the size of its ionic model
1199  // and plugins, and the associated displacments between the nodes
1200  SF::vector<int> impsize(numNode), impdsp;
1201  for(int i=0; i<numNode; i++)
1202  impsize[i] = imp_mem[(int)IIFmask[i]];
1203  SF::dsp_from_cnt(impsize, impdsp);
1204 
1205  // this also gives us the size of the total buffer
1206  size_t num_entr = SF::sum(impsize);
1207  SF::vector<char> impdata(num_entr);
1208 
1209  // now we can memcpy the ionic models from the IonIfBase arrays into our buffer
1210  for (int imp_idx = 0; imp_idx < this->N_IIF; imp_idx++)
1211  {
1212  IonIfBase* iif = this->IIF[imp_idx];
1213  // When using data layout optimization, checkpointing indices have to be
1214  // change a bit. Without DLO, this value is simply 1
1215  std::size_t vec_size = iif->get_type().dlo_vector_size();
1216  for (int imp_nod_idx = 0; imp_nod_idx < this->N_Nodes[imp_idx]; imp_nod_idx+=vec_size) {
1217  int index = imp_nod_idx / vec_size;
1218  int loc = this->NodeLists[imp_idx][imp_nod_idx];
1219  int svSize = iif->get_sv_size();
1220 
1221  char* write = impdata.data() + impdsp[loc];
1222  char* read = (char*) iif->get_sv_address() + index * svSize;
1223  memcpy(write, read, svSize);
1224 
1225  int shift = svSize;
1226  for (auto& plug : iif->plugins()) {
1227  int plugsize = plug->get_sv_size();
1228 
1229  write = impdata.data() + impdsp[loc] + shift;
1230  read = (char*) plug->get_sv_address() + index * plugsize;
1231  memcpy(write, read, plugsize);
1232  shift += plugsize;
1233  }
1234  }
1235  }
1236 
1237  // finally we can write out the IonIf data to disk
1238  SF::root_write_ordered<int, char>(fd, canord, impsize.data(), impdata.data(),
1239  numNode, num_entr, PETSC_COMM_WORLD);
1240  delete [] canord;
1241  delete [] imp_mem;
1242 
1243  if(fd) fclose(fd);
1244 
1245  double dump_time = timing(t1, t0);
1246  log_msg(logger, 0, 0, " in %.3f seconds.\n", dump_time);
1247 }
1248 
1271 float MULTI_IF::restore_state(const char *fname, mesh_t gid, bool close)
1272 {
1273  static FILE_SPEC in = NULL;
1274  static char *last_fn = NULL;
1275  double t0, t1;
1276 
1277  t0 = get_time();
1278  int error = 0;
1279  int my_rank = get_rank();
1280  int mpi_size = get_size();
1281 
1282  if (my_rank == 0) {
1283  if (fname) in = f_open(fname, "r");
1284 
1285  if (!in) {
1286  if (fname) log_msg(logger, 5, 0, "Error: cannot open file: %s\n", fname);
1287  else log_msg(logger, 5, 0, "Error: file stream not open yet");
1288  error++;
1289  }
1290  else {
1291  if (strcmp(read_bin_string(in), Magic_MIIF_ID) ) {
1292  log_msg(logger, 5, 0, "%s is not a recognized MIIF dump file", fname);
1293  error++;
1294  }
1295  }
1296  }
1297 
1298  if (get_global(error, MPI_SUM))
1299  EXIT(1);
1300 
1301  if (fname) {
1302  free(last_fn);
1303  last_fn = strdup(fname);
1304  }
1305  else {
1306  fname = last_fn;
1307  }
1308 
1309  unsigned int format, version;
1310  float time = 0.0f;
1311  time_t save_date;
1312  f_read_par(&format, sizeof(unsigned int), 1, in);
1313  f_read_par(&version, sizeof(unsigned int), 1, in);
1314  f_read_par(&save_date, sizeof(time_t), 1, in);
1315  f_read_par(&time, sizeof(float), 1, in);
1316  log_msg(logger, 0, 0, "Restoring time %f from %s (format v%d) generated\n\tby calling "
1317  "program r%d on %s", time, fname, format, version, ctime(&save_date) );
1318 
1319  int savedNum, glob_numNode;
1320  MPI_Allreduce(&this->numNode, &glob_numNode, 1, MPI_INT, MPI_SUM, PETSC_COMM_WORLD);
1321 
1322  f_read_par(&savedNum, sizeof(int), 1, in);
1323  if (savedNum != glob_numNode) {
1324  log_msg(logger, 5, 0, "expecting %d nodes but read %d nodes", glob_numNode, savedNum);
1325  EXIT(1);
1326  }
1327 
1328  int num_gdata;
1329  f_read_par(&num_gdata, sizeof(int), 1, in);
1330  SF::scattering & petsc2canon = *get_permutation(gid, PETSC_TO_CANONICAL, 1);
1331 
1332  for (int g = 0; g < num_gdata; g++) {
1333  char *datatype = read_bin_string_par(in);
1334  int i;
1335  for (i = 0; i < NUM_IMP_DATA_TYPES; i++) {
1336  if (!strcmp(datatype, imp_data_names[i]) ) {
1337  if (this->gdata[i]) {
1338  log_msg(logger, 0, 0, "\tRestoring global data %s", imp_data_names[i], datatype);
1339 
1340  FILE* fd = my_rank == 0 ? in->fd : NULL;
1341  this->gdata[i]->read_binary<SF_real>(fd);
1342  bool fwd = false;
1343  petsc2canon(*this->gdata[i], fwd);
1344  } else {
1345  log_msg(logger, 0, 0, "\tGlobal data %s not used", datatype);
1346  }
1347  break;
1348  }
1349  }
1350  if (i == NUM_IMP_DATA_TYPES)
1351  log_msg(logger, 3, 0, "\tSaved global data %s not recognized", datatype);
1352 
1353  if ( (my_rank == 0) && ((i == NUM_IMP_DATA_TYPES) || !this->gdata[i]) )
1354  fseek(in->fd, this->gdata[i]->gsize() * sizeof(SF_real), SEEK_CUR);
1355 
1356  free(datatype);
1357  }
1358 
1359  int N_IIF;
1360  f_read_par(&N_IIF, sizeof(int), 1, in);
1361 
1362  // read in the IMPs for each region
1363  IMPinfo *imps = static_cast<IMPinfo *>(calloc(N_IIF, sizeof(IMPinfo)));
1364  int *IMPsz = static_cast<int *>(malloc(N_IIF * sizeof(int)));
1365 
1366  for (int i = 0; i < N_IIF; i++) {
1367  imps[i].name = read_bin_string_par(in);
1368  f_read_par(&imps[i].sz, sizeof(int), 1, in);
1369  f_read_par(&imps[i].nplug, sizeof(int), 1, in);
1370 
1371  imps[i].offset = 0;
1372  IMPsz[i] = imps[i].sz;
1373  imps[i].plug = static_cast<IMPinfo *>(calloc(sizeof(IMPinfo), imps[i].nplug));
1374 
1375  for (int j = 0; j < imps[i].nplug; j++) {
1376  imps[i].plug[j].offset = IMPsz[i];
1377  imps[i].plug[j].name = read_bin_string_par(in);
1378  f_read_par(&imps[i].plug[j].sz, sizeof(int), 1, in);
1379  IMPsz[i] += imps[i].plug[j].sz;
1380  }
1381  }
1382 
1383  // compare new and old regions to determine which can be copied
1384  for (int i = 0; i < N_IIF; i++) {
1385  if (i >= this->N_IIF) {
1386  log_msg(NULL, 3, 0, "Saved IMP region %d out of range", i);
1387  continue;
1388  }
1389  if (strcmp(imps[i].name, this->IIF[i]->get_type().get_name().c_str()) ) {
1390  log_msg(NULL, 3, 0, "Saved IMP region %d ionic model does not match that of IMP region %d", i, i);
1391  continue;
1392  }
1393  if (imps[i].sz != this->IIF[i]->get_sv_size()) {
1394  log_msg(NULL, 3, 0, "Saved IMP region %d size does not match current IMP size", i);
1395  continue;
1396  }
1397 
1398  // the IM is good to go
1399  imps[i].compatible = true;
1400  for (int j = 0; j < imps[i].nplug; j++)
1401  for (int k = 0; k < this->IIF[i]->plugins().size(); k++)
1402  if (!strcmp(imps[i].plug[j].name, this->IIF[i]->plugins()[k]->get_type().get_name().c_str()) &&
1403  (imps[i].plug[j].sz == this->IIF[i]->plugins()[k]->get_sv_size()) ) {
1404  log_msg(NULL, 3, 0, "Saved IMP region %d plugin %s compatible", i, imps[i].plug[j].name);
1405  imps[i].plug[j].map = k;
1406  imps[i].plug[j].compatible = true;
1407  break;
1408  }
1409 
1410  }
1411 
1412  // read in region mask
1413  IIF_Mask_t *canMask = static_cast<char *>(malloc(glob_numNode*sizeof(this->IIFmask[0]) ));
1414  f_read_par(canMask, sizeof(IIF_Mask_t), glob_numNode, in);
1415 
1416  // determine relative IMP data offsets based on canonical ordering
1417  size_t *offset = static_cast<size_t *>(calloc(glob_numNode+1, sizeof(size_t)));
1418  for (int i = 1; i <= glob_numNode; i++)
1419  offset[i] = offset[i-1] + IMPsz[static_cast<int>(canMask[i-1])];
1420 
1421  // read in the state variable data
1422  // TODO: Aurel: Ideally rank0 would read in the data and communicate it, but right now I dont have
1423  // time to code this up. Thus, at least the ranks will read their data one by one so
1424  // we dont saturate the file system.
1425  long SVstart;
1426  if (my_rank == 0) {
1427  SVstart = ftell(in->fd);
1428  f_close(in);
1429  }
1430 
1431  MPI_Bcast(&SVstart, sizeof(long), MPI_BYTE, 0, PETSC_COMM_WORLD);
1432  int mismatch = 0;
1433 
1434  // we need to map from a local algebraic index to a global canonical index
1435  const sf_mesh & mesh = get_mesh(gid);
1436  const SF::vector<mesh_int_t> & canon_nbr = mesh.get_numbering(SF::NBR_SUBMESH);
1437  const SF::vector<mesh_int_t> & alg_idx = mesh.pl.algebraic_nodes();
1438  SF::vector<mesh_int_t> loc2canon(alg_idx.size());
1439 
1440  for(size_t i=0; i<alg_idx.size(); i++) loc2canon[i] = canon_nbr[alg_idx[i]];
1441 
1442  // the ranks read the file one-by-one
1443  for (int pid = 0; pid < mpi_size; pid++) {
1444  if (my_rank == pid) {
1445  if(in) delete in;
1446  in = f_open(fname, "r");
1447 
1448  for (int i = 0; i < N_IIF; i++) {
1449  fseek(in->fd, SVstart, SEEK_SET);
1450  mismatch += this->IIF[i]->restore(in, this->N_Nodes[i], this->NodeLists[i],
1451  canMask, offset, imps+i, loc2canon.data());
1452  }
1453 
1454  f_close(in);
1455  }
1456  MPI_Barrier(PETSC_COMM_WORLD);
1457  }
1458  free(canMask);
1459  free(IMPsz);
1460 
1461  MPI_Reduce(my_rank ? &mismatch : MPI_IN_PLACE, &mismatch, 1, MPI_INT, MPI_SUM, 0, PETSC_COMM_WORLD);
1462  if ( (my_rank == 0) && mismatch)
1463  log_msg(NULL, 3, 0, "Number of nonmatching nodes: %d", mismatch);
1464 
1465  if ( (my_rank == 0) && (!close) ) {
1466  if(in) delete in;
1467  in = f_open(fname, "r");
1468 
1469  // position at end of current MIIF data
1470  fseek(in->fd, SVstart+offset[glob_numNode], SEEK_SET);
1471  }
1472 
1473  free(offset);
1474  for (int i = 0; i < N_IIF; i++)
1475  free(imps[i].plug);
1476  free(imps);
1477 
1478  double restore_time = timing(t1, t0);
1479  log_msg(logger, 0, 0, "State restored from file %s in %.3f seconds.\n", fname, restore_time);
1480 
1481  return time;
1482 }
1483 
1495 void MULTI_IF::sv_dump_add_by_name_list( int region, char *imp_name,
1496  char *reg_name, char *sv_lst, char *plg_lst,
1497  char *plg_sv_lst, double t, double dump_dt) {
1498  char file[8000], svs[1024], plgs[1024], plgsvs[1024];
1499  char *e, *l, *p, *svnames, *plgnames, *plgsvnames;
1500 
1501  strcpy(svs, sv_lst);
1502  strcpy(plgs, plg_lst);
1503  strcpy(plgsvs, plg_sv_lst);
1504  svnames = &svs[0];
1505  plgnames = &plgs[0];
1506  plgsvnames = &plgsvs[0];
1507 
1508  // override default setting for dump interval
1509  this->svd.intv = dump_dt;
1510 
1511  // override default start time for dumping (in case of a restart)
1512  this->svd.t_dump = t;
1513 
1514  // parse ionic model sv list
1515  while ( (e = get_next_list(svnames, ',') ) ) {
1516  snprintf(file, sizeof file, "%s.%s.bin", reg_name, svnames);
1517  this->sv_dump_add_by_name(region, imp_name, svnames, reg_name, file);
1518  svnames = e;
1519  }
1520 
1521  // Information is provided as follows:
1522  // region[X].plugin = "PLG_A:PLG_B"
1523  // region[X].plug_sv_dumps = "m,h,n:x,y"
1524 
1525  // parse plugin list
1526  while ( (p = get_next_list(plgnames, ':') ) ) {
1527  // get sv list corresponding to current plugin
1528  l = get_next_list(plgsvnames, ':');
1529  while ( (e = get_next_list(plgsvnames, ',') ) ) {
1530  snprintf(file, sizeof file, "%s_%s.%s.bin", reg_name, plgnames, plgsvnames);
1531  this->sv_dump_add_by_name(region, plgnames, plgsvnames, reg_name, file);
1532  plgsvnames = e;
1533  }
1534  plgnames = p;
1535  plgsvnames = l;
1536  }
1537 } // sv_dump_add_by_name_list
1538 
1554 void MULTI_IF::sv_dump_add(int region, const IonType& type, int offset, int size, int dtype,
1555  const char *filename, const char *regname) {
1556  IonIfBase *IF = this->IIF[region];
1557  int n = this->svd.n;
1558 
1559  // find the proper IMP
1560  if (IF->get_type() != type) {
1561  int i;
1562  for (i = 0; i < IF->plugins().size(); i++)
1563  if (IF->plugins()[i]->get_type() == type)
1564  break;
1565  if (i == IF->plugins().size()) {
1566  log_msg(logger, 2, 0, "Warning: IMP %s not found in Region %s\n",
1567  type.get_name().c_str(), regname);
1568  return;
1569  } else {
1570  IF = IF->plugins()[i];
1571  }
1572  }
1573 
1574  this->svd.active = 1;
1575  this->svd.n++;
1576  this->svd.hdls = static_cast<FILE_SPEC *>(realloc(this->svd.hdls, this->svd.n*sizeof(this->svd.hdls[0])));
1577  this->svd.fn = static_cast<char **>(realloc(this->svd.fn, this->svd.n*sizeof(char *)));
1578  this->svd.reg = static_cast<int *>(realloc(this->svd.reg, this->svd.n*sizeof(int)));
1579  this->svd.svnames = static_cast<char **>(realloc(this->svd.svnames, this->svd.n*sizeof(char *)));
1580  this->svd.n_dumps = 0;
1581  this->svd.nwr = 0;
1582  this->svd.offset = static_cast<int *>(realloc(this->svd.offset, this->svd.n*sizeof(int)));
1583  this->svd.size = static_cast<int *>(realloc(this->svd.size, this->svd.n*sizeof(int)));
1584  this->svd.dtype = static_cast<int *>(realloc(this->svd.dtype, this->svd.n*sizeof(int)));
1585  this->svd.svtab = static_cast<void **>(realloc(this->svd.svtab, this->svd.n*sizeof(void *)));
1586  this->svd.svsize = static_cast<size_t *>(realloc(this->svd.svsize, this->svd.n*sizeof(size_t)));
1587  this->svd.num = static_cast<int *>(realloc(this->svd.num, this->svd.n*sizeof(int)));
1588  this->svd.reg = static_cast<int *>(realloc(this->svd.reg, this->svd.n*sizeof(int)));
1589 
1590  if (!get_rank() ) {
1591 #ifdef HDF5
1592  assert(0);
1593 #else // ifdef HDF5
1594  this->svd.hdls[n] = f_open(filename, "w+");
1595 #endif // ifdef HDF5
1596  } else {
1597  this->svd.hdls[n] = NULL;
1598  }
1599 
1600  this->svd.fn[n] = dupstr(filename);
1601  this->svd.reg[n] = region;
1602  this->svd.svnames[n] = NULL; // has to be filled in in calling routine
1603  this->svd.offset[n] = offset;
1604  this->svd.size[n] = size;
1605  this->svd.dtype[n] = dtype;
1606  this->svd.svtab[n] = IF->get_sv_address();
1607  this->svd.svsize[n] = IF->get_sv_size();
1608  this->svd.num[n] = IF->get_num_node();
1609 } // sv_dump_add
1610 
1621 int MULTI_IF::sv_dump_add_by_name(int region, char *impname,
1622  char *svname, char *regname, char *filename) {
1623  int offset, size, dtype, added = 0;
1624  char *svtypename = NULL;
1625  char *filename_bin;
1626 
1627  IonType* type = get_ion_type(std::string(impname));
1628  filename_bin = strcat(filename, ".bin");
1629 
1630  assert(type != NULL);
1631  if (type->get_sv_offset(svname, &offset, &size) ) {
1632  type->get_sv_type(svname, &dtype, &svtypename);
1633  this->sv_dump_add(region, *type, offset, size, dtype, filename_bin, regname);
1634  this->svd.svnames[this->svd.n-1] = dupstr(svname);
1635  added = 1;
1636  } else {
1637  log_msg(NULL, 1, 0, "No state variable added to dump list");
1638  }
1639 
1640  return added;
1641 }
1642 
1655 {
1656  assert(miif->gdata[Vm] != NULL);
1657 
1658  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++) {
1659  for (int n = 0; n < miif->N_IIF; n++)
1660  if (USED_DAT(miif->IIF[n], imp_data_flag[i]) && (miif->gdata[i] == NULL)) {
1661  SF::init_vector(&miif->gdata[i], miif->gdata[Vm]);
1662  break;
1663  }
1664  }
1665 }
1666 
1675  int num_mech_data = 0;
1676 
1677  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++) {
1678  if ( (Lambda_DATA_FLAG == i) || (delLambda_DATA_FLAG == i) ||
1679  (Tension_DATA_FLAG == i) || (tension_component_DATA_FLAG == i) ) {
1680  for (int n = 0; n < this->N_IIF; n++)
1681  if (USED_DAT(this->IIF[n], imp_data_flag[i]) )
1682  num_mech_data++;
1683  }
1684  }
1685  return static_cast<bool>(num_mech_data);
1686 }
1687 
1705 void MULTI_IF::transmem_stim_species(float curr, const char *species,
1706  float beta, int *node, int numnode) {
1707  float z;
1708  if (strstr(species, "Ca") != NULL) {
1709  curr *= 1000; // Ca uses micromolar instead of millimolar
1710  z = 2;
1711  } else if (strstr(species, "Na") != NULL) {
1712  z = 1;
1713  } else if (strstr(species, "Cl") != NULL) {
1714  z = -1;
1715  } else if (strstr(species, "K") != NULL) {
1716  z = 1;
1717  } else {
1718  log_msg(logger, 5, 0, "Unimplemented ion species: %s\n", species);
1719  exit(1);
1720  }
1721  int offset, sz;
1722  SVgetfcn ion_get;
1723  SVputfcn ion_put;
1724  static int warned = 0;
1725 
1726  int rank = get_rank();
1727  SF::vector<int> layout;
1728  layout_from_count(numnode, layout, PETSC_COMM_WORLD);
1729  int my_low_idx = layout[rank], my_high_idx = layout[rank+1];
1730 
1731  for (int i = 0; i < numnode; i++) {
1732  if ((node[i] < my_low_idx) || (node[i] >= my_high_idx)) continue; // not on processor
1733 
1734  for (int j = 0; j < this->N_IIF; j++) {
1735  if (this->NodeLists[j] == NULL) continue; // no stim nodes on processor
1736 
1737  if ((node[i] < this->NodeLists[j][0]) ||
1738  (node[i] > this->NodeLists[j][this->N_Nodes[j]-1]) ) {
1739  continue; // not in this IMP
1740  } else {
1741  ion_get = this->IIF[j]->get_type().get_sv_offset(species, &offset, &sz);
1742  if (ion_get == NULL) {
1743  if (!warned) {
1744  warned = 1;
1745  log_msg(logger, 2, 0, "Ion species not present in ionic model: %s\n",
1746  species);
1747  }
1748  return;
1749  }
1750  ion_put = getPutSV(ion_get);
1751  }
1752 
1753  // does this model provide a correct conversion factor?
1754  double stim_curr;
1755  cell_geom g = this->IIF[j]->cgeom();
1756  if (g.sl_i2c != NDEF) {
1757  stim_curr = curr*g.sl_i2c/z;
1758  } else {
1759  // conversion with generic factor
1760  stim_curr = curr*beta/FARADAY/z*1.e-2; // convert to millimolar/L
1761  }
1762 
1763  for (int k = 0; k < this->N_Nodes[j]; k++)
1764  if (node[i] == this->NodeLists[j][k])
1765  ion_put(*this->IIF[j], k, offset, ion_get(*this->IIF[j], k, offset)-stim_curr);
1766  }
1767  }
1768 } // transmem_stim_species
1769 
1776 void MULTI_IF::MIIF_change_dt( double Dt) {
1777  this->dt = Dt;
1778  for (int i = 0; i < this->N_IIF; i++) {
1779  this->IIF[i]->destroy_luts();
1780  this->IIF[i]->set_dt((float) Dt);
1781  this->iontypes[i].get().construct_tables(*this->IIF[i]);
1782 
1783  for (int j = 0; j < this->IIF[i]->plugins().size(); j++) {
1784  auto& plugin = this->IIF[i]->plugins()[j];
1785  plugin->destroy_luts();
1786  plugin->set_dt((float) Dt);
1787  this->plugtypes[i][j].get().construct_tables(*plugin);
1788  }
1789  }
1790 }
1791 
1792 #define MEMFREE(A) free(A)
1793 
1801  assert(m->doppel);
1802 
1803  for (int i = 0; i < m->N_IIF; i++) {
1804  m->IIF[i]->get_type().destroy_ion_if(m->IIF[i]);
1805  }
1806 
1807  for (int j = 0; j < NUM_IMP_DATA_TYPES; j++) {
1808  if (m->gdata[j])
1809  delete m->gdata[j];
1810  }
1811 
1812  m->zero_data();
1813 }
1814 
1815 #undef MEMFREE
1816 
1824 void doppel_update(MULTI_IF *orig, MULTI_IF *miif_doppel) {
1825  // removed assertion, need to be able to update bidirectionally
1826  // assert(doppel->doppel);
1827 
1828  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++)
1829  if (orig->gdata[i])
1830  SF::init_vector(&miif_doppel->gdata[i], orig->gdata[i]);
1831 
1832  for (int i = 0; i < orig->N_IIF; i++) {
1833  miif_doppel->IIF[i]->copy_SVs_from(*orig->IIF[i], false);
1834 
1835  for (int j = 0; j < orig->IIF[i]->plugins().size(); j++)
1836  miif_doppel->IIF[i]->plugins()[j]->copy_SVs_from(*orig->IIF[i]->plugins()[j], false);
1837  }
1838 
1839  miif_doppel->getRealData();
1840 }
1841 
1852 void doppel_MIIF(MULTI_IF *orig, MULTI_IF *miif_doppel) {
1853  *miif_doppel = *orig;
1854  miif_doppel->doppel = true;
1855 
1856  // deal with the local portion of the global data
1857  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++)
1858  if (orig->gdata[i])
1859  SF::init_vector(&miif_doppel->gdata[i], orig->gdata[i]);
1860 
1861  // copy the state variables
1862  miif_doppel->IIF = {};
1863 
1864  for (int i = 0; i < orig->N_IIF; i++) {
1865  // Make a new IonIf for this clone and copy state variables and plugins
1866  miif_doppel->IIF.push_back(orig->IIF[i]->get_type().make_ion_if(orig->IIF[i]->get_target(),
1867  orig->IIF[i]->get_num_node(), orig->plugtypes[i]));
1868  miif_doppel->IIF[i]->copy_SVs_from(*orig->IIF[i], true);
1869  miif_doppel->IIF[i]->copy_plugins_from(*orig->IIF[i]);
1870  }
1871 }
1872 
1880 bool isIMPdata(const char *sv) {
1881  return IMPdataLabel2Index(sv) != -1;
1882 }
1883 
1891 int IMPdataLabel2Index(const char *sv) {
1892  int imp_data_id = -1;
1893 
1894  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++) {
1895  if (strcmp(sv, imp_data_names[i]) == 0) {
1896  imp_data_id = i;
1897  break;
1898  }
1899  }
1900  return imp_data_id;
1901 }
1902 
1903 /* copy all of the IMP data from one node to another
1904  *
1905  * \param IF ionic model
1906  * \param from local IMP number
1907  * \param to local IMP number
1908  */
1909 void dup_IMP_node_state(IonIfBase& IF, int from, int to, GlobalData_t **localdata) {
1910  for (int i = 0; i < NUM_IMP_DATA_TYPES; i++)
1911  if (localdata && localdata[i])
1912  localdata[i][to] = localdata[i][from];
1913 
1914  IF.for_each([&](IonIfBase& imp) {
1915  char **sv_list;
1916  int sv_list_size = imp.get_type().get_sv_list(&sv_list);
1917 
1918  for (int i = 0; i < sv_list_size; i++) {
1919  char *sv_name = sv_list[i];
1920  int sv_offset;
1921  int sv_size;
1922  SVgetfcn sv_get = imp.get_type().get_sv_offset(sv_name, &sv_offset, &sv_size);
1923  if (sv_get == NULL) {
1924  throw std::runtime_error(std::string(__func__) + " error: " + sv_name + " is not a valid state variable for the " + imp.get_type().get_name() + " model.");
1925  }
1926  SVputfcn sv_put = getPutSV(sv_get);
1927 
1928  GlobalData_t sv_val = sv_get(imp, from, sv_offset);
1929  sv_put(imp, to, sv_offset, sv_val);
1930  }
1931  free(sv_list);
1932  });
1933 }
1934 
1935 } // namespace limpet
#define NDEF
definition of cell geometry
Definition: ION_IF.h:112
#define FILENAME_BUF
#define MAX_TRACE_LINE_LEN
Define multiple ionic models to be used in different regions.
#define FARADAY
Faraday's constant.
Definition: MULTI_ION_IF.h:147
#define USED_DAT(I, F)
Definition: MULTI_ION_IF.h:300
int mesh_int_t
Definition: SF_container.h:46
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
#define SLIST_APPEND(S, P)
Definition: basics.h:52
#define FLUSH
Definition: basics.h:311
#define LOCAL
Definition: basics.h:309
#define fmemopen_
Definition: basics.h:433
#define STRUCT_ZERO(S)
Definition: basics.h:51
virtual S * ptr()=0
virtual void release_ptr(S *&p)=0
size_t read_binary(FILE *fd)
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false)=0
virtual void add_scaled(const abstract_vector< T, S > &vec, S k)=0
size_t write_binary(FILE *fd)
Write a vector to HD in binary. File descriptor is already set up.
bool in_b(const T idx)
return whether idx is in set B
Definition: SF_container.h:355
T backward_map(T idx) const
Map one index from b to a.
Definition: SF_container.h:263
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:419
Container for a PETSc VecScatter.
void forward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Forward scattering.
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
T * data()
Pointer to the vector's start.
Definition: SF_vector.h:91
Represents the ionic model and plug-in (IMP) data structure.
Definition: ION_IF.h:138
const IonType & get_type() const
Gets this IMP's model type.
Definition: ION_IF.cc:138
std::vector< IonIfBase * > & plugins()
Returns a vector containing the plugins of this IMP.
Definition: ION_IF.cc:178
ts & get_tstp()
Gets the time stepper.
Definition: ION_IF.cc:202
void for_each(const std::function< void(IonIfBase &)> &consumer)
Executes the consumer functions on this IMP and each of its plugins.
Definition: ION_IF.cc:411
IonIfBase * parent() const
Gets the parent IMP.
Definition: ION_IF.cc:170
int get_num_node() const
Gets the number of nodes handled by this IMP.
Definition: ION_IF.cc:142
int miifIdx
imp index within miif
Definition: ION_IF.h:150
virtual std::size_t get_sv_size() const =0
Gets the size of the structure this IMP uses for state variables.
virtual void * get_sv_address()=0
Gets the raw address of the state variables for this IMP.
void compute(int start, int end, GlobalData_t **data)
Perform ionic model computation for 1 time step.
Definition: ION_IF.cc:264
Abstract class representing an ionic model type.
Definition: ion_type.h:59
bool is_plugin() const
Returns whether this model is a plugin or not.
Definition: ion_type.cc:27
const std::string & get_name() const
Gets the model name.
Definition: ion_type.cc:23
virtual SVgetfcn get_sv_offset(const char *svname, int *off, int *sz) const =0
Get the offset and size of a state variable of the model, as well as an access function.
virtual int get_sv_type(const char *svname, int *type, char **type_name) const =0
Determines the type of a SV.
virtual size_t dlo_vector_size() const =0
Gets the vector size when using data layout optimization (DLO).
int numNode
local number of nodes
Definition: MULTI_ION_IF.h:209
bool extUpdateVm
flag indicating update function for Vm
Definition: MULTI_ION_IF.h:207
int dump_svs(opencarp::base_timer *)
std::vector< IonIfBase * > IIF
array of IIF's
Definition: MULTI_ION_IF.h:201
opencarp::sf_vec * gdata[NUM_IMP_DATA_TYPES]
data used by all IMPs
Definition: MULTI_ION_IF.h:215
void sv_dump_add_by_name_list(int, char *, char *, char *, char *, char *, double, double)
std::vector< Target > targets
target for each region
Definition: MULTI_ION_IF.h:212
std::vector< IonTypeList > plugtypes
plugins types for each region
Definition: MULTI_ION_IF.h:214
IonTypeList iontypes
type for each region
Definition: MULTI_ION_IF.h:211
void dump_state(char *, float, opencarp::mesh_t gid, bool, unsigned int)
SV_DUMP svd
state variable dump
Definition: MULTI_ION_IF.h:202
void sv_dump_add(int, const IonType &, int, int, int, const char *, const char *)
void transmem_stim_species(float, const char *, float, int *, int)
void initialize_currents(double, int)
GlobalData_t *** ldata
data local to each IMP
Definition: MULTI_ION_IF.h:204
int adjust_MIIF_variables(const char *variable, const SF::vector< int > &indices, const SF::vector< double > &values)
float restore_state(const char *, opencarp::mesh_t gid, bool)
int N_IIF
how many different IIF's
Definition: MULTI_ION_IF.h:210
void compute_ionic_current(bool flag_send=1, bool flag_receive=1)
GPU kernel to emulate the add_scaled call made to adjust the Vm values when the update to Vm is not m...
int sv_dump_add_by_name(int, char *, char *, char *, char *)
Trace_Info * trace_info
Information about traces.
Definition: MULTI_ION_IF.h:206
double dt
time step (ms)
Definition: MULTI_ION_IF.h:218
int * N_Nodes
#nodes for each IMP
Definition: MULTI_ION_IF.h:199
int ** NodeLists
local partitioned node lists for each IMP stored
Definition: MULTI_ION_IF.h:200
GlobalData_t * procdata[NUM_IMP_DATA_TYPES]
data for this processor
Definition: MULTI_ION_IF.h:203
void dump_luts_MIIF(bool)
bool * contiguous
whether a region is contiguously numbered
Definition: MULTI_ION_IF.h:205
opencarp::FILE_SPEC logger
Definition: MULTI_ION_IF.h:216
int numSubDt
number of sub-dt time steps
Definition: MULTI_ION_IF.h:219
void MIIF_change_dt(double)
IIF_Mask_t * IIFmask
region for each node
Definition: MULTI_ION_IF.h:213
std::string name
name for MIIF region
Definition: MULTI_ION_IF.h:197
bool doppel
is this a shallow clone?
Definition: MULTI_ION_IF.h:198
void releaseRealDataDuringInit()
void dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
Definition: SF_vector.h:310
T sum(const vector< T > &vec)
Compute sum of a vector's entries.
Definition: SF_vector.h:340
void local_petsc_to_nodal_mapping(const meshdata< T, S > &mesh, index_mapping< T > &petsc_to_nodal)
void layout_from_count(const T count, vector< T > &layout, MPI_Comm comm)
Definition: SF_network.h:201
void init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:122
@ NBR_SUBMESH
Submesh nodal numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:191
void doppel_MIIF(MULTI_IF *orig, MULTI_IF *miif_doppel)
int getGlobalNodalIndex(IonIfBase &pIF, int relIdx)
void initializeIMPData(MULTI_IF *pMIIF)
int determine_write_ranges(int N, size_t *offset, size_t bufsize, int **ranges)
int get_plug_flag(char *plgstr, int *out_num_plugins, IonTypeList &out_plugins)
GlobalData_t(* SVgetfcn)(IonIfBase &, int, int)
Definition: ion_type.h:48
SVputfcn getPutSV(SVgetfcn)
void doppel_update(MULTI_IF *orig, MULTI_IF *miif_doppel)
double Real
Definition: MULTI_ION_IF.h:151
void dup_IMP_node_state(IonIfBase &IF, int from, int to, GlobalData_t **localdata)
int should_print_bounds_exceeded_messages()
@ CPU
baseline CPU model generated with the original opencarp code generator
Definition: target.h:48
IonType * get_ion_type(const std::string &name)
int IMPdataLabel2Index(const char *sv)
void CreateIIFGlobalNodeLsts(MULTI_IF *pMIIF)
SF_real GlobalData_t
Definition: limpet_types.h:27
void allocate_shared_data(MULTI_IF *)
bool isIMPdata(const char *)
bool is_gpu(Target const target)
Checks if this is a GPU target.
Definition: target.cc:64
int current_global_node(int local_node)
void initialize_params_MIIF(MULTI_IF *pMIIF)
void(* SVputfcn)(IonIfBase &, int, int, GlobalData_t)
Definition: ion_type.h:49
char * get_sv(void *tab, int offset, int n, int svSize, int size)
void alloc_MIIF(MULTI_IF *pMIIF)
void freeIMPData(MULTI_IF *pMIIF)
std::vector< std::reference_wrapper< IonType > > IonTypeList
Definition: ion_type.h:288
constexpr T max(T a, T b)
Definition: ion_type.h:31
void deallocate_on_target(Target target, T *ptr)
Utility function for deallocating memory on a target. See TargetAllocator.
Definition: target.h:314
void CreateIIFNodeLsts_(int, char *, int **, int ***, int)
void CreateIIFLocalNodeLsts(MULTI_IF *pMIIF)
void initialize_ionic_IF(MULTI_IF *pMIIF)
void update_ts(ts *ptstp)
Definition: ION_IF.cc:460
float current_global_time()
void close_trace(MULTI_IF *MIIF)
void dump_trace(MULTI_IF *MIIF, limpet::Real time)
char * get_next_list(char *lst, char delimiter)
Definition: ION_IF.cc:506
int int_cmp(const void *a, const void *b)
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.
char IIF_Mask_t
Definition: ion_type.h:50
void free_doppel(MULTI_IF *m)
char * tokstr_r(char *s1, const char *s2, char **lasts)
Definition: ION_IF.cc:82
int set_print_bounds_exceeded_messages(int newval)
SF::abstract_vector< mesh_int_t, double > sf_vec
Definition: sf_interface.h:49
sf_mesh & get_mesh(const mesh_t gt)
Get a mesh by specifying the gridID.
Definition: sf_interface.cc:33
SF::scattering * get_permutation(const int mesh_id, const int perm_id, const int dpn)
Get the PETSC to canonical permutation scattering for a given mesh and number of dpn.
SF::meshdata< mesh_int_t, mesh_real_t > sf_mesh
Definition: sf_interface.h:47
char * read_bin_string(FILE_SPEC in)
Definition: basics.cc:228
void f_read_par(void *ptr, size_t size, size_t nmemb, FILE_SPEC stream, MPI_Comm comm)
Parallel fread. Root reads, then broadcasts.
Definition: basics.cc:171
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:276
T get_global(T in, MPI_Op OP, MPI_Comm comm=PETSC_COMM_WORLD)
Do a global reduction on a variable.
Definition: basics.h:233
void write_bin_string(FILE_SPEC out, const char *s)
Definition: basics.cc:219
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:135
char * read_bin_string_par(FILE_SPEC in)
Definition: basics.cc:239
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:60
void get_time(double &tm)
Definition: basics.h:436
int get_size(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:290
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
file_desc * FILE_SPEC
Definition: basics.h:138
#define PETSC_TO_CANONICAL
Permute algebraic data from PETSC to canonical ordering.
Definition: sf_interface.h:76
int offset
offset into node data
Definition: ION_IF.h:132
int sz
storage required
Definition: ION_IF.h:127
char * name
IMP name.
Definition: ION_IF.h:126
int map
which plugin does this IMO match
Definition: ION_IF.h:131
IMPinfo * plug
plugins
Definition: ION_IF.h:129
bool compatible
does IM match stored IM
Definition: ION_IF.h:130
int nplug
number of plugins
Definition: ION_IF.h:128
data structure to manage state variable file dumps
Definition: MULTI_ION_IF.h:160
size_t * svsize
state variable sizes
Definition: MULTI_ION_IF.h:176
char ** fn
array to store file names
Definition: MULTI_ION_IF.h:164
double intv
time interval for sv dumps
Definition: MULTI_ION_IF.h:167
int * dtype
data type
Definition: MULTI_ION_IF.h:173
int n_dumps
keep track of number of dumped time slices
Definition: MULTI_ION_IF.h:169
long nwr
keep track of number of written tokens
Definition: MULTI_ION_IF.h:170
int * size
sizes of SV to dump
Definition: MULTI_ION_IF.h:172
void ** svtab
state variable tables
Definition: MULTI_ION_IF.h:175
int * reg
array to store region ids
Definition: MULTI_ION_IF.h:165
double t_dump
next instant for sv dump
Definition: MULTI_ION_IF.h:168
int * num
number of nodes
Definition: MULTI_ION_IF.h:174
char ** svnames
array to store sv names
Definition: MULTI_ION_IF.h:166
int n
#state variables we want to dump
Definition: MULTI_ION_IF.h:162
int * offset
offsets into structure for SV
Definition: MULTI_ION_IF.h:171
opencarp::FILE_SPEC * hdls
array of file handles to sv output files
Definition: MULTI_ION_IF.h:163
data structure to manage trace dumps. Should eventually be combined with the state variable dumps,...
Definition: MULTI_ION_IF.h:183
bool ignored
globally not found
Definition: MULTI_ION_IF.h:185
bool found
found on this node
Definition: MULTI_ION_IF.h:184
opencarp::FILE_SPEC file
Definition: MULTI_ION_IF.h:188
int node_idx
local node number
Definition: MULTI_ION_IF.h:186
float sl_i2c
convert sl-currents in uA/cm^2 to mM/L without valence
Definition: ION_IF.h:118
int cnt
Definition: ION_IF.h:86