NEURON
nrn_setup.cpp
Go to the documentation of this file.
1 /*
2 # =============================================================================
3 # Copyright (c) 2016 - 2021 Blue Brain Project/EPFL
4 #
5 # See top-level LICENSE file for details.
6 # =============================================================================.
7 */
8 
9 #include <algorithm>
10 #include <vector>
11 #include <map>
12 #include <cstring>
13 #include <mutex>
14 
16 #include "coreneuron/nrnconf.h"
26 #include "coreneuron/mpi/nrnmpi.h"
35 #include "coreneuron/io/phase1.hpp"
36 #include "coreneuron/io/phase2.hpp"
37 #include "coreneuron/io/phase3.hpp"
40 
41 // callbacks into nrn/src/nrniv/nrnbbcore_write.cpp
44 
45 
46 /// --> Coreneuron
50 
51 void (*nrn2core_group_ids_)(int*);
52 
53 extern "C" {
54 SetupTransferInfo* (*nrn2core_get_partrans_setup_info_)(int ngroup,
55  int cn_nthread,
56  size_t cn_sidt_size);
57 }
58 
60  int& bsize,
61  int& n_pr,
62  void**& vpr,
63  int& n_trajec,
64  int*& types,
65  int*& indices,
66  double**& pvars,
67  double**& varrays);
68 
69 void (*nrn2core_trajectory_values_)(int tid, int n_pr, void** vpr, double t);
70 
71 void (*nrn2core_trajectory_return_)(int tid, int n_pr, int bsize, int vecsz, void** vpr, double t);
72 
73 int (*nrn2core_all_spike_vectors_return_)(std::vector<double>& spikevec, std::vector<int>& gidvec);
74 
75 void (*nrn2core_all_weights_return_)(std::vector<double*>& weights);
76 
77 // file format defined in cooperation with nrncore/src/nrniv/nrnbbcore_write.cpp
78 // single integers are ascii one per line. arrays are binary int or double
79 // Note that regardless of the gid contents of a group, since all gids are
80 // globally unique, a filename convention which involves the first gid
81 // from the group is adequate. Also note that balance is carried out from a
82 // per group perspective and launching a process consists of specifying
83 // a list of group ids (first gid of the group) for each process.
84 //
85 // <firstgid>_1.dat
86 // n_presyn, n_netcon
87 // output_gids (npresyn) with -(type+1000*index) for those acell with no gid
88 // netcon_srcgid (nnetcon) -(type+1000*index) refers to acell with no gid
89 // -1 means the netcon has no source (not implemented)
90 // Note that the negative gids are only thread unique and not process unique.
91 // We create a thread specific hash table for the negative gids for each thread
92 // when <firstgid>_1.dat is read and then destroy it after <firstgid>_2.dat
93 // is finished using it. An earlier implementation which attempted to
94 // encode the thread number into the negative gid
95 // (i.e -ith - nth*(type +1000*index)) failed due to not large enough
96 // integer domain size.
97 // Note that for file transfer it is an error if a negative srcgid is
98 // not in the same thread as the target. This is because there it may
99 // not be the case that threads in a NEURON process end up on same process
100 // in CoreNEURON. NEURON will raise an error if this
101 // is the case. However, for direct memory transfer, it is allowed that
102 // a negative srcgid may be in a different thread than the target. So
103 // nrn2core_get_dat1 has a last arg netcon_negsrcgid_tid that specifies
104 // for the negative gids in netcon_srcgid (in that order) the source thread.
105 //
106 // <firstgid>_2.dat
107 // n_real_cell, n_output, n_real_output, nnode
108 // ndiam - 0 if no mechanism has dparam with diam semantics, or nnode
109 // nmech - includes artcell mechanisms
110 // for the nmech tml mechanisms
111 // type, nodecount
112 // nidata, nvdata, nweight
113 // v_parent_index (nnode)
114 // actual_a, b, area, v (nnode)
115 // diam - if ndiam > 0. Note that only valid diam is for those nodes with diam semantics mechanisms
116 // for the nmech tml mechanisms
117 // nodeindices (nodecount) but only if not an artificial cell
118 // data (nodecount*param_size)
119 // pdata (nodecount*dparam_size) but only if dparam_size > 0 on this side.
120 // output_vindex (n_presyn) >= 0 associated with voltages -(type+1000*index) for acell
121 // output_threshold (n_real_output)
122 // netcon_pnttype (nnetcon) <=0 if a NetCon does not have a target.
123 // netcon_pntindex (nnetcon)
124 // weights (nweight)
125 // delays (nnetcon)
126 // for the nmech tml mechanisms that have a nrn_bbcore_write method
127 // type
128 // icnt
129 // dcnt
130 // int array (number specified by the nodecount nrn_bbcore_write
131 // to be intepreted by this side's nrn_bbcore_read method)
132 // double array
133 // #VectorPlay_instances, for each of these instances
134 // 4 (VecPlayContinuousType)
135 // mtype
136 // index (from Memb_list.data)
137 // vecsize
138 // yvec
139 // tvec
140 //
141 // The critical issue requiring careful attention is that a coreneuron
142 // process reads many coreneuron thread files with a result that, although
143 // the conceptual
144 // total n_pre is the sum of all the n_presyn from each thread as is the
145 // total number of output_gid, the number of InputPreSyn instances must
146 // be computed here from a knowledge of all thread's netcon_srcgid after
147 // all thread's output_gids have been registered. We want to save the
148 // "individual allocation of many small objects" memory overhead by
149 // allocating a single InputPreSyn array for the entire process.
150 // For this reason cellgroup data are divided into two separate
151 // files with the first containing output_gids and netcon_srcgid which are
152 // stored in the nt.presyns array and nt.netcons array respectively
153 namespace coreneuron {
154 static OMP_Mutex mut;
155 
156 /// Vector of maps for negative presyns
157 std::vector<std::map<int, PreSyn*>> neg_gid2out;
158 /// Maps for ouput and input presyns
159 std::map<int, PreSyn*> gid2out;
160 std::map<int, InputPreSyn*> gid2in;
161 
162 /// InputPreSyn.nc_index_ to + InputPreSyn.nc_cnt_ give the NetCon*
163 std::vector<NetCon*> netcon_in_presyn_order_;
164 
165 /// Only for setup vector of netcon source gids
166 std::vector<int*> nrnthreads_netcon_srcgid;
167 
168 /// If a nrnthreads_netcon_srcgid is negative, need to determine the thread when
169 /// in order to use the correct neg_gid2out[tid] map
170 std::vector<std::vector<int>> nrnthreads_netcon_negsrcgid_tid;
171 
172 /* read files.dat file and distribute cellgroups to all mpi ranks */
173 void nrn_read_filesdat(int& ngrp, int*& grp, const char* filesdat) {
174  patstimtype = nrn_get_mechtype("PatternStim");
177  grp = new int[ngrp + 1];
178  (*nrn2core_group_ids_)(grp);
179  return;
180  }
181 
182  FILE* fp = fopen(filesdat, "r");
183 
184  if (!fp) {
185  nrn_fatal_error("No input file ( %s ) with nrnthreads, exiting...", filesdat);
186  }
187 
188  char version[256];
189  nrn_assert(fscanf(fp, "%s\n", version) == 1);
191 
192  int iNumFiles = 0;
193  nrn_assert(fscanf(fp, "%d\n", &iNumFiles) == 1);
194 
195  // temporary strategem to figure out if model uses gap junctions while
196  // being backward compatible
197  if (iNumFiles == -1) {
198  nrn_assert(fscanf(fp, "%d\n", &iNumFiles) == 1);
199  nrn_have_gaps = true;
200  if (nrnmpi_myid == 0) {
201  printf("Model uses gap junctions\n");
202  }
203  }
204 
205  if (nrnmpi_numprocs > iNumFiles && nrnmpi_myid == 0) {
206  printf(
207  "Info : The number of input datasets are less than ranks, some ranks will be idle!\n");
208  }
209 
210  ngrp = 0;
211  grp = new int[iNumFiles / nrnmpi_numprocs + 1];
212 
213  // irerate over gids in files.dat
214  for (int iNum = 0; iNum < iNumFiles; ++iNum) {
215  int iFile;
216 
217  nrn_assert(fscanf(fp, "%d\n", &iFile) == 1);
218  if ((iNum % nrnmpi_numprocs) == nrnmpi_myid) {
219  // A "-1" entry means that this rank should not be assigned further gid groups.
220  // It is a way to create files.dat files which deterministically assign gid groups to
221  // ranks, particularly useful for very large simulations which required load balancing.
222  if (iFile == -1) {
223  break;
224  }
225  grp[ngrp] = iFile;
226  ngrp++;
227  }
228  }
229 
230  fclose(fp);
231 }
232 
233 void netpar_tid_gid2ps(int tid, int gid, PreSyn** ps, InputPreSyn** psi) {
234  /// for gid < 0 returns the PreSyn* in the thread (tid) specific map.
235  *ps = nullptr;
236  *psi = nullptr;
237 
238  if (gid >= 0) {
239  auto gid2out_it = gid2out.find(gid);
240  if (gid2out_it != gid2out.end()) {
241  *ps = gid2out_it->second;
242  } else {
243  auto gid2in_it = gid2in.find(gid);
244  if (gid2in_it != gid2in.end()) {
245  *psi = gid2in_it->second;
246  }
247  }
248  } else {
249  auto gid2out_it = neg_gid2out[tid].find(gid);
250  if (gid2out_it != neg_gid2out[tid].end()) {
251  *ps = gid2out_it->second;
252  }
253  }
254 }
255 
257  // allocate the process wide InputPreSyn array
258  // all the output_gid have been registered and associated with PreSyn.
259  // now count the needed InputPreSyn by filling the netpar::gid2in map
260  gid2in.clear();
261 
262  // now have to fill the new table
263  // do not need to worry about negative gid overlap since only use
264  // it to search for PreSyn in this thread.
265 
266  std::vector<InputPreSyn*> inputpresyn_;
267 
268  for (int ith = 0; ith < nrn_nthread; ++ith) {
269  NrnThread& nt = nrn_threads[ith];
270  // associate gid with InputPreSyn and increase PreSyn and InputPreSyn count
271  nt.n_input_presyn = 0;
272  // if single thread or file transfer then definitely empty.
273  std::vector<int>& negsrcgid_tid = nrnthreads_netcon_negsrcgid_tid[ith];
274  size_t i_tid = 0;
275  for (int i = 0; i < nt.n_netcon; ++i) {
276  int gid = nrnthreads_netcon_srcgid[ith][i];
277  if (gid >= 0) {
278  /// If PreSyn or InputPreSyn is already in the map
279  auto gid2out_it = gid2out.find(gid);
280  if (gid2out_it != gid2out.end()) {
281  /// Increase PreSyn count
282  ++gid2out_it->second->nc_cnt_;
283  continue;
284  }
285  auto gid2in_it = gid2in.find(gid);
286  if (gid2in_it != gid2in.end()) {
287  /// Increase InputPreSyn count
288  ++gid2in_it->second->nc_cnt_;
289  continue;
290  }
291 
292  /// Create InputPreSyn and increase its count
293  InputPreSyn* psi = new InputPreSyn;
294  ++psi->nc_cnt_;
295  gid2in[gid] = psi;
296  inputpresyn_.push_back(psi);
297  ++nt.n_input_presyn;
298  } else {
299  int tid = nt.id;
300  if (!negsrcgid_tid.empty()) {
301  tid = negsrcgid_tid[i_tid++];
302  }
303  auto gid2out_it = neg_gid2out[tid].find(gid);
304  if (gid2out_it != neg_gid2out[tid].end()) {
305  /// Increase negative PreSyn count
306  ++gid2out_it->second->nc_cnt_;
307  }
308  }
309  }
310  }
311 
312  // now, we can opportunistically create the NetCon* pointer array
313  // to save some memory overhead for
314  // "large number of small array allocation" by
315  // counting the number of NetCons each PreSyn and InputPreSyn point to.
316  // Conceivably the nt.netcons could become a process global array
317  // in which case the NetCon* pointer array could become an integer index
318  // array. More speculatively, the index array could be eliminated itself
319  // if the process global NetCon array were ordered properly but that
320  // would interleave NetCon from different threads. Not a problem for
321  // serial threads but the reordering would propagate to nt.pntprocs
322  // if the NetCon data pointers are also replaced by integer indices.
323 
324  // First, allocate the pointer array.
325  int n_nc = 0;
326  for (int ith = 0; ith < nrn_nthread; ++ith) {
327  n_nc += nrn_threads[ith].n_netcon;
328  }
329  netcon_in_presyn_order_.resize(n_nc);
330  n_nc = 0;
331 
332  // fill the indices with the offset values and reset the nc_cnt_
333  // such that we use the nc_cnt_ in the following loop to assign the NetCon
334  // to the right place
335  // for PreSyn
336  int offset = 0;
337  for (int ith = 0; ith < nrn_nthread; ++ith) {
338  NrnThread& nt = nrn_threads[ith];
339  for (int i = 0; i < nt.n_presyn; ++i) {
340  PreSyn& ps = nt.presyns[i];
341  ps.nc_index_ = offset;
342  offset += ps.nc_cnt_;
343  ps.nc_cnt_ = 0;
344  }
345  }
346  // for InputPreSyn
347  for (auto psi: inputpresyn_) {
348  psi->nc_index_ = offset;
349  offset += psi->nc_cnt_;
350  psi->nc_cnt_ = 0;
351  }
352 
353  inputpresyn_.clear();
354 
355  // with gid to InputPreSyn and PreSyn maps we can setup the multisend
356  // target lists.
357  if (use_multisend_) {
358 #if NRN_MULTISEND
360 #endif
361  }
362 
363  // fill the netcon_in_presyn_order and recompute nc_cnt_
364  // note that not all netcon_in_presyn will be filled if there are netcon
365  // with no presyn (ie. nrnthreads_netcon_srcgid[nt.id][i] = -1) but that is ok since they are
366  // only used via ps.nc_index_ and ps.nc_cnt_;
367  for (int ith = 0; ith < nrn_nthread; ++ith) {
368  NrnThread& nt = nrn_threads[ith];
369  // if single thread or file transfer then definitely empty.
370  std::vector<int>& negsrcgid_tid = nrnthreads_netcon_negsrcgid_tid[ith];
371  size_t i_tid = 0;
372  for (int i = 0; i < nt.n_netcon; ++i) {
373  NetCon* nc = nt.netcons + i;
374  int gid = nrnthreads_netcon_srcgid[ith][i];
375  int tid = ith;
376  if (!negsrcgid_tid.empty() && gid < -1) {
377  tid = negsrcgid_tid[i_tid++];
378  }
379  PreSyn* ps;
380  InputPreSyn* psi;
381  netpar_tid_gid2ps(tid, gid, &ps, &psi);
382  if (ps) {
383  netcon_in_presyn_order_[ps->nc_index_ + ps->nc_cnt_] = nc;
384  ++ps->nc_cnt_;
385  ++n_nc;
386  } else if (psi) {
387  netcon_in_presyn_order_[psi->nc_index_ + psi->nc_cnt_] = nc;
388  ++psi->nc_cnt_;
389  ++n_nc;
390  }
391  }
392  }
393 
394  /// Resize the vector to its actual size of the netcons put in it
395  netcon_in_presyn_order_.resize(n_nc);
396 }
397 
398 /// Clean up
400  for (int ith = 0; ith < nrn_nthread; ++ith) {
401  if (nrnthreads_netcon_srcgid[ith])
402  delete[] nrnthreads_netcon_srcgid[ith];
403  }
404  nrnthreads_netcon_srcgid.clear();
406  neg_gid2out.clear();
407 }
408 
409 void nrn_setup(const char* filesdat,
410  bool is_mapping_needed,
411  CheckPoints& checkPoints,
412  bool run_setup_cleanup,
413  const char* datpath,
414  const char* restore_path,
415  double* mindelay) {
416  double time = nrn_wtime();
417 
418  int ngroup;
419  int* gidgroups;
420  nrn_read_filesdat(ngroup, gidgroups, filesdat);
421  UserParams userParams(ngroup,
422  gidgroups,
423  datpath,
424  strlen(restore_path) == 0 ? datpath : restore_path,
425  checkPoints);
426 
427 
428  // temporary bug work around. If any process has multiple threads, no
429  // process can have a single thread. So, for now, if one thread, make two.
430  // Fortunately, empty threads work fine.
431  // Allocate NrnThread* nrn_threads of size ngroup (minimum 2)
432  // Note that rank with 0 dataset/cellgroup works fine
433  nrn_threads_create(userParams.ngroup <= 1 ? 2 : userParams.ngroup);
434 
435  // from nrn_has_net_event create pnttype2presyn for use in phase2.
436  auto& memb_func = corenrn.get_memb_funcs();
437  auto& pnttype2presyn = corenrn.get_pnttype2presyn();
439  pnttype2presyn.clear();
440  pnttype2presyn.resize(memb_func.size(), -1);
441  for (size_t i = 0; i < nrn_has_net_event_.size(); ++i) {
442  pnttype2presyn[nrn_has_net_event_[i]] = i;
443  }
444 
446 
447  if (nrn_nthread > 1) {
448  // NetCvode construction assumed one thread. Need nrn_nthread instances
449  // of NetCvodeThreadData. Here since possible checkpoint restore of
450  // tqueue at end of phase2.
451  nrn_p_construct();
452  }
453 
454  if (use_solve_interleave) {
456  }
457 
458  /// Reserve vector of maps of size ngroup for negative gid-s
459  /// std::vector< std::map<int, PreSyn*> > neg_gid2out;
460  neg_gid2out.resize(userParams.ngroup);
461 
462  // bug fix. gid2out is cumulative over all threads and so do not
463  // know how many there are til after phase1
464  // A process's complete set of output gids and allocation of each thread's
465  // nt.presyns and nt.netcons arrays.
466  // Generates the gid2out map which is needed
467  // to later count the required number of InputPreSyn
468  /// gid2out - map of output presyn-s
469  /// std::map<int, PreSyn*> gid2out;
470  gid2out.clear();
471 
473  for (int i = 0; i < nrn_nthread; ++i)
474  nrnthreads_netcon_srcgid[i] = nullptr;
475 
476  // Gap junctions used to be done first in the sense of reading files
477  // and calling gap_mpi_setup. But during phase2, gap_thread_setup and
478  // gap_indices_permute were called after NrnThread.data was in its final
479  // layout and mechanism permutation was determined. This is no longer
480  // ideal as it necessitates keeping setup_info_ in existence to the end
481  // of phase2. So gap junction setup is deferred to after phase2.
482 
484  if (corenrn_file_mode) {
485  coreneuron::phase_wrapper<coreneuron::phase::one>(userParams, !corenrn_file_mode);
486  } else {
488  Phase1 p1{n->id};
489  p1.populate(*n, mut);
490  });
491  }
492 
493  // from the gid2out map and the nrnthreads_netcon_srcgid array,
494  // fill the gid2in, and from the number of entries,
495  // allocate the process wide InputPreSyn array
497 
498  // read the rest of the gidgroup's data and complete the setup for each
499  // thread.
500  /* nrn_multithread_job supports serial, pthread, and openmp. */
501  coreneuron::phase_wrapper<coreneuron::phase::two>(userParams, !corenrn_file_mode);
502 
503  // gap junctions
504  // Gaps are done after phase2, in order to use layout and permutation
505  // information via calls to legacy_index2pointer.
506  if (nrn_have_gaps) {
508  if (!corenrn_embedded) {
510  coreneuron::phase_wrapper<coreneuron::gap>(userParams);
511  } else {
512  nrn_partrans::setup_info_ = (*nrn2core_get_partrans_setup_info_)(userParams.ngroup,
513  nrn_nthread,
514  sizeof(sgid_t));
515  }
516 
519 
520  // Whether allocated in NEURON or here, delete here.
521  delete[] nrn_partrans::setup_info_;
522  nrn_partrans::setup_info_ = nullptr;
523  }
524 
525  if (is_mapping_needed)
526  coreneuron::phase_wrapper<coreneuron::phase::three>(userParams, !corenrn_file_mode);
527 
528  *mindelay = set_mindelay(*mindelay);
529 
530  if (run_setup_cleanup) // if run_setup_cleanup==false, user must call nrn_setup_cleanup() later
532 
533 #if INTERLEAVE_DEBUG
534  // mk_cell_indices debug code is supposed to be used with cell-per-core permutations
536  mk_cell_indices();
537  }
538 #endif
539 
540  /// Allocate memory for fast_imem calculation
542 
543  /// Generally, tables depend on a few parameters. And if those parameters change,
544  /// then the table needs to be recomputed. This is obviously important in NEURON
545  /// since the user can change those parameters at any time. However, there is no
546  /// c example for CoreNEURON so can't see what it looks like in that context.
547  /// Boils down to setting up a function pointer of the function _check_table_thread(),
548  /// which is only executed by StochKV.c.
549  nrn_mk_table_check(); // was done in nrn_thread_memblist_setup in multicore.c
550 
551  size_t model_size_bytes;
552 
555  model_size_bytes = model_size(true);
556  } else {
557  model_size_bytes = model_size(false);
558  }
559 
560  if (nrnmpi_myid == 0 && !corenrn_param.is_quiet()) {
561  printf(" Setup Done : %.2lf seconds \n", nrn_wtime() - time);
562 
563  if (model_size_bytes < 1024) {
564  printf(" Model size : %ld bytes\n", model_size_bytes);
565  } else if (model_size_bytes < 1024 * 1024) {
566  printf(" Model size : %.2lf kB\n", model_size_bytes / 1024.);
567  } else if (model_size_bytes < 1024 * 1024 * 1024) {
568  printf(" Model size : %.2lf MB\n", model_size_bytes / (1024. * 1024.));
569  } else {
570  printf(" Model size : %.2lf GB\n", model_size_bytes / (1024. * 1024. * 1024.));
571  }
572  }
573 
574  delete[] userParams.gidgroups;
575 }
576 
578  for (NrnThreadMembList* tml = nt.tml; tml; tml = tml->next) {
579  Memb_func& mf = corenrn.get_memb_func(tml->index);
580  Memb_list* ml = tml->ml;
581  if (mf.thread_size_) {
583  if (mf.thread_mem_init_) {
584  {
585  const std::lock_guard<OMP_Mutex> lock(mut);
586  (*mf.thread_mem_init_)(ml->_thread);
587  }
588  }
589  } else {
590  ml->_thread = nullptr;
591  }
592  }
593 }
594 
595 void read_phasegap(NrnThread& nt, UserParams& userParams) {
596  auto& F = userParams.file_reader[nt.id];
597  if (F.fail()) {
598  return;
599  }
600 
601  F.checkpoint(0);
602 
603  int sidt_size = F.read_int();
604  assert(sidt_size == int(sizeof(sgid_t)));
605  std::size_t ntar = F.read_int();
606  std::size_t nsrc = F.read_int();
607 
608  auto& si = nrn_partrans::setup_info_[nt.id];
609  si.src_sid.resize(nsrc);
610  si.src_type.resize(nsrc);
611  si.src_index.resize(nsrc);
612  if (nsrc) {
613  F.read_array<sgid_t>(si.src_sid.data(), nsrc);
614  F.read_array<int>(si.src_type.data(), nsrc);
615  F.read_array<int>(si.src_index.data(), nsrc);
616  }
617 
618  si.tar_sid.resize(ntar);
619  si.tar_type.resize(ntar);
620  si.tar_index.resize(ntar);
621  if (ntar) {
622  F.read_array<sgid_t>(si.tar_sid.data(), ntar);
623  F.read_array<int>(si.tar_type.data(), ntar);
624  F.read_array<int>(si.tar_index.data(), ntar);
625  }
626 
627 #if CORENRN_DEBUG
628  printf("%d read_phasegap tid=%d nsrc=%d ntar=%d\n", nrnmpi_myid, nt.id, nsrc, ntar);
629  for (int i = 0; i < nsrc; ++i) {
630  printf("src %z %d %d\n", size_t(si.src_sid[i]), si.src_type[i], si.src_index[i]);
631  }
632  for (int i = 0; i < ntar; ++i) {
633  printf("tar %z %d %d\n", size_t(si.src_sid[i]), si.src_type[i], si.src_index[i]);
634  }
635 #endif
636 }
637 
638 // This function is related to nrn_dblpntr2nrncore in Neuron to determine which values should
639 // be transferred from CoreNeuron. Types correspond to the value to be transferred based on
640 // mech_type enum or non-artificial cell mechanisms.
641 // take into account alignment, layout, permutation
642 // only voltage, i_membrane_ or mechanism data index allowed. (mtype 0 means time)
643 double* legacy_index2pointer(int mtype, int index, NrnThread& nt) {
644  if (mtype == voltage) { // voltage
645  int ix = index; // relative to _actual_v
646  nrn_assert((ix >= 0) && (ix < nt.end));
647  if (nt._permute) {
648  node_permute(&ix, 1, nt._permute);
649  }
650  return nt._actual_v + ix;
651  } else if (mtype == i_membrane_) { // membrane current from fast_imem calculation
652  int ix = index; // relative to nrn_fast_imem->nrn_sav_rhs
653  nrn_assert((ix >= 0) && (ix < nt.end));
654  if (nt._permute) {
655  node_permute(&ix, 1, nt._permute);
656  }
657  return nt.nrn_fast_imem->nrn_sav_rhs + ix;
658  } else if (mtype > 0 && mtype < static_cast<int>(corenrn.get_memb_funcs().size())) { //
659  Memb_list* ml = nt._ml_list[mtype];
660  nrn_assert(ml);
661 
662  const std::vector<int>& array_dims = corenrn.get_array_dims()[mtype];
663  int padded_node_count = nrn_soa_padded_size(ml->nodecount, Layout::SoA);
664 
665  auto soaos_index = legacy2soaos_index(index, array_dims);
666  auto cnrn_index =
667  soaos2cnrn_index(soaos_index, array_dims, padded_node_count, ml->_permute);
668 
669  return ml->data + cnrn_index;
670  } else if (mtype == 0) { // time
671  return &nt._t;
672  } else {
673  printf("legacy_index2pointer does not handle mtype=%d\n", mtype);
674  nrn_assert(0);
675  }
676  return nullptr;
677 }
678 
679 // from i to (icnt, isz)
680 void nrn_inverse_i_layout(int i, int& icnt, int cnt, int& isz, int sz, int layout) {
681  if (layout == Layout::AoS) {
682  icnt = i / sz;
683  isz = i % sz;
684  } else if (layout == Layout::SoA) {
685  int padded_cnt = nrn_soa_padded_size(cnt, layout);
686  icnt = i % padded_cnt;
687  isz = i / padded_cnt;
688  } else {
689  assert(0);
690  }
691 }
692 
693 /**
694  * Cleanup global ion map created during mechanism registration
695  *
696  * In case of coreneuron standalone execution nrn_ion_global_map
697  * can be deleted at the end of execution. But in case embedded
698  * run via neuron, mechanisms are registered only once i.e. during
699  * first call to coreneuron. This is why we call cleanup only in
700  * case of standalone coreneuron execution via nrniv-core or
701  * special-core.
702  *
703  * @todo coreneuron should have finalise callback which can be
704  * called from NEURON for final memory cleanup including global
705  * state like registered mechanisms and ions map.
706  */
708  for (int i = 0; i < nrn_ion_global_map_size; i++) {
710  }
712  nrn_ion_global_map = nullptr;
714 }
715 
717  delete[] std::exchange(nt._fornetcon_perm_indices, nullptr);
718  delete[] std::exchange(nt._fornetcon_weight_perm, nullptr);
719 }
720 
721 /* nrn_threads_free() presumes all NrnThread and NrnThreadMembList data is
722  * allocated with malloc(). This is not the case here, so let's try and fix
723  * things up first. */
724 
725 void nrn_cleanup() {
726  clear_event_queue(); // delete left-over TQItem
727  for (auto psi: gid2in) {
728  delete psi.second;
729  }
730  gid2in.clear();
731  gid2out.clear();
732 
733  // clean nrnthread_chkpnt
734  if (nrnthread_chkpnt) {
735  delete[] nrnthread_chkpnt;
736  nrnthread_chkpnt = nullptr;
737  }
738 
739  // clean NrnThreads
740  for (int it = 0; it < nrn_nthread; ++it) {
741  NrnThread* nt = nrn_threads + it;
742  NrnThreadMembList* next_tml = nullptr;
745  for (NrnThreadMembList* tml = nt->tml; tml; tml = next_tml) {
746  Memb_list* ml = tml->ml;
747 
748  mod_f_t s = corenrn.get_memb_func(tml->index).destructor;
749  if (s) {
750  (*s)(nt, ml, tml->index);
751  }
752 
753  // Moved from below as priv_dtor is now deleting the RANDOM streams,
754  // and at this moment need an undeleted pdata.
755  // Destroy the global variables struct allocated in nrn_init
756  if (auto* const priv_dtor = corenrn.get_memb_func(tml->index).private_destructor) {
757  (*priv_dtor)(nt, ml, tml->index);
758  assert(!ml->instance);
759  assert(!ml->global_variables);
760  assert(ml->global_variables_size == 0);
761  }
762 
763  ml->data = nullptr; // this was pointing into memory owned by nt
764  free_memory(ml->pdata);
765  ml->pdata = nullptr;
767  ml->nodeindices = nullptr;
768  if (ml->_permute) {
769  delete[] ml->_permute;
770  ml->_permute = nullptr;
771  }
772 
773  if (ml->_thread) {
774  free_memory(ml->_thread);
775  ml->_thread = nullptr;
776  }
777 
779  if (nrb) {
780  if (nrb->_size) {
781  free_memory(nrb->_pnt_index);
783  free_memory(nrb->_nrb_t);
784  free_memory(nrb->_nrb_flag);
785  free_memory(nrb->_displ);
786  free_memory(nrb->_nrb_index);
787  }
788  free_memory(nrb);
789  ml->_net_receive_buffer = nullptr;
790  }
791 
793  if (nsb) {
794  delete nsb;
795  ml->_net_send_buffer = nullptr;
796  }
797 
798  if (tml->dependencies)
799  free(tml->dependencies);
800 
801  next_tml = tml->next;
802  free_memory(tml->ml);
803  free_memory(tml);
804  }
805 
806  nt->_actual_rhs = nullptr;
807  nt->_actual_d = nullptr;
808  nt->_actual_a = nullptr;
809  nt->_actual_b = nullptr;
810 
812  nt->_v_parent_index = nullptr;
813 
814  free_memory(nt->_data);
815  nt->_data = nullptr;
816 
817  free(nt->_idata);
818  nt->_idata = nullptr;
819 
820  free_memory(nt->_vdata);
821  nt->_vdata = nullptr;
822 
823  if (nt->_permute) {
824  delete[] nt->_permute;
825  nt->_permute = nullptr;
826  }
827 
828  if (nt->presyns_helper) {
830  nt->presyns_helper = nullptr;
831  }
832 
833  delete[] std::exchange(nt->pntprocs, nullptr);
834  delete[] std::exchange(nt->presyns, nullptr);
835 
836  if (nt->pnt2presyn_ix) {
837  for (size_t i = 0; i < corenrn.get_has_net_event().size(); ++i) {
838  if (nt->pnt2presyn_ix[i]) {
839  free(nt->pnt2presyn_ix[i]);
840  }
841  }
843  }
844 
845  if (nt->netcons) {
846  delete[] nt->netcons;
847  nt->netcons = nullptr;
848  }
849 
850  if (nt->weights) {
851  free_memory(nt->weights);
852  nt->weights = nullptr;
853  }
854 
855  if (nt->_shadow_rhs) {
857  nt->_shadow_rhs = nullptr;
858  }
859 
860  if (nt->_shadow_d) {
861  free_memory(nt->_shadow_d);
862  nt->_shadow_d = nullptr;
863  }
864 
865  if (nt->_net_send_buffer_size) {
867  nt->_net_send_buffer = nullptr;
868  nt->_net_send_buffer_size = 0;
869  }
870 
871  if (nt->_watch_types) {
872  free(nt->_watch_types);
873  nt->_watch_types = nullptr;
874  }
875 
876  // mapping information is available only for non-empty NrnThread
877  if (nt->mapping && nt->ncell) {
878  delete ((NrnThreadMappingInfo*) nt->mapping);
879  }
880 
881  free_memory(nt->_ml_list);
882 
883  if (nt->nrn_fast_imem) {
884  fast_imem_free();
885  }
886  }
887 
888 #if NRN_MULTISEND
890 #endif
891 
892  netcon_in_presyn_order_.clear();
893 
895 
896  if (!corenrn.get_pnttype2presyn().empty()) {
897  corenrn.get_pnttype2presyn().clear();
898  }
899 
901 
903 }
904 
906  if (nt.trajec_requests) {
908  if (tr->n_trajec) {
909  delete[] tr->vpr;
910  if (tr->scatter) {
911  delete[] tr->scatter;
912  }
913  if (tr->varrays) {
914  delete[] tr->varrays;
915  }
916  delete[] tr->gather;
917  }
918  delete nt.trajec_requests;
919  nt.trajec_requests = nullptr;
920  }
921 }
922 
923 void read_phase1(NrnThread& nt, UserParams& userParams) {
924  Phase1 p1{userParams.file_reader[nt.id]};
925 
926  // Protect gid2in, gid2out and neg_gid2out
927  p1.populate(nt, mut);
928 }
929 
930 void read_phase2(NrnThread& nt, UserParams& userParams) {
931  Phase2 p2;
933  p2.read_direct(nt.id, nt);
934  } else {
935  p2.read_file(userParams.file_reader[nt.id], nt);
936  }
937  p2.populate(nt, userParams);
938 }
939 
940 /** read mapping information for neurons */
941 void read_phase3(NrnThread& nt, UserParams& userParams) {
942  /** mapping information for all neurons in single NrnThread */
943  NrnThreadMappingInfo* ntmapping = new NrnThreadMappingInfo();
944 
945  Phase3 p3;
947  p3.read_direct(ntmapping);
948  } else {
949  auto& F = userParams.file_reader[nt.id];
950  F.restore_checkpoint();
951  p3.read_file(F, ntmapping);
952  }
953 
954  // make number #cells match with mapping size
955  nrn_assert((int) ntmapping->size() == nt.ncell);
956 
957  // set pointer in NrnThread
958  nt.mapping = (void*) ntmapping;
959  nt.summation_report_handler_ = std::make_unique<SummationReportMapping>();
960 }
961 
962 /* Returns the size of the dynamically allocated memory for NrnThreadMembList
963  * Includes:
964  * - Size of NrnThreadMembList
965  * - Size of Memb_list
966  * - Size of nodeindices
967  * - Size of _permute
968  * - Size of _thread
969  * - Size of NetReceive and NetSend Buffers
970  * - Size of int variables
971  * - Size of double variables (If include_data is enabled. Those variables are already counted
972  * since they point to nt->_data.)
973  */
974 size_t memb_list_size(NrnThreadMembList* tml, bool include_data) {
975  size_t nbyte = sizeof(NrnThreadMembList) + sizeof(Memb_list);
976  nbyte += tml->ml->nodecount * sizeof(int);
977  if (tml->ml->_permute) {
978  nbyte += tml->ml->nodecount * sizeof(int);
979  }
980  if (tml->ml->_thread) {
981  Memb_func& mf = corenrn.get_memb_func(tml->index);
982  nbyte += mf.thread_size_ * sizeof(ThreadDatum);
983  }
984  if (tml->ml->_net_receive_buffer) {
985  nbyte += sizeof(NetReceiveBuffer_t) + tml->ml->_net_receive_buffer->size_of_object();
986  }
987  if (tml->ml->_net_send_buffer) {
988  nbyte += sizeof(NetSendBuffer_t) + tml->ml->_net_send_buffer->size_of_object();
989  }
990  if (include_data) {
991  nbyte += corenrn.get_prop_param_size()[tml->index] * tml->ml->nodecount * sizeof(double);
992  }
993  nbyte += corenrn.get_prop_dparam_size()[tml->index] * tml->ml->nodecount * sizeof(Datum);
994 #ifdef DEBUG
995  int i = tml->index;
996  printf("%s %d psize=%d ppsize=%d cnt=%d nbyte=%ld\n",
997  corenrn.get_memb_func(i).sym,
998  i,
1001  tml->ml->nodecount,
1002  nbyte);
1003 #endif
1004  return nbyte;
1005 }
1006 
1007 /// Approximate count of number of bytes for the gid2out map
1008 size_t output_presyn_size(void) {
1009  if (gid2out.empty()) {
1010  return 0;
1011  }
1012  size_t nbyte = sizeof(gid2out) + sizeof(int) * gid2out.size() +
1013  sizeof(PreSyn*) * gid2out.size();
1014 #ifdef DEBUG
1015  printf(" gid2out table bytes=~%ld size=%ld\n", nbyte, gid2out.size());
1016 #endif
1017  return nbyte;
1018 }
1019 
1020 size_t input_presyn_size(void) {
1021  if (gid2in.empty()) {
1022  return 0;
1023  }
1024  size_t nbyte = sizeof(gid2in) + sizeof(int) * gid2in.size() +
1025  sizeof(InputPreSyn*) * gid2in.size();
1026 #ifdef DEBUG
1027  printf(" gid2in table bytes=~%ld size=%ld\n", nbyte, gid2in.size());
1028 #endif
1029  return nbyte;
1030 }
1031 
1032 size_t model_size(bool detailed_report) {
1033  long nbyte = 0;
1034  size_t sz_nrnThread = sizeof(NrnThread);
1035  size_t sz_presyn = sizeof(PreSyn);
1036  size_t sz_input_presyn = sizeof(InputPreSyn);
1037  size_t sz_netcon = sizeof(NetCon);
1038  size_t sz_pntproc = sizeof(Point_process);
1039  size_t nccnt = 0;
1040 
1041  std::vector<long> size_data(13, 0);
1042  std::vector<long> global_size_data_min(13, 0);
1043  std::vector<long> global_size_data_max(13, 0);
1044  std::vector<long> global_size_data_sum(13, 0);
1045  std::vector<float> global_size_data_avg(13, 0.0);
1046 
1047  for (int i = 0; i < nrn_nthread; ++i) {
1048  NrnThread& nt = nrn_threads[i];
1049  size_t nb_nt = 0; // per thread
1050  nccnt += nt.n_netcon;
1051 
1052  // Memb_list size
1053  int nmech = 0;
1054  for (auto tml = nt.tml; tml; tml = tml->next) {
1055  nb_nt += memb_list_size(tml, false);
1056  ++nmech;
1057  }
1058 
1059  // basic thread size includes mechanism data and G*V=I matrix
1060  nb_nt += sz_nrnThread;
1061  nb_nt += nt._ndata * sizeof(double) + nt._nidata * sizeof(int) + nt._nvdata * sizeof(void*);
1062  nb_nt += nt.end * sizeof(int); // _v_parent_index
1063 
1064  // network connectivity
1065  nb_nt += nt.n_pntproc * sz_pntproc + nt.n_netcon * sz_netcon + nt.n_presyn * sz_presyn +
1066  nt.n_input_presyn * sz_input_presyn + nt.n_weight * sizeof(double);
1067  nbyte += nb_nt;
1068 
1069 #ifdef DEBUG
1070  printf("ncell=%d end=%d nmech=%d\n", nt.ncell, nt.end, nmech);
1071  printf("ndata=%ld nidata=%ld nvdata=%ld\n", nt._ndata, nt._nidata, nt._nvdata);
1072  printf("nbyte so far %ld\n", nb_nt);
1073  printf("n_presyn = %d sz=%ld nbyte=%ld\n", nt.n_presyn, sz_presyn, nt.n_presyn * sz_presyn);
1074  printf("n_input_presyn = %d sz=%ld nbyte=%ld\n",
1075  nt.n_input_presyn,
1076  sz_input_presyn,
1077  nt.n_input_presyn * sz_input_presyn);
1078  printf("n_pntproc=%d sz=%ld nbyte=%ld\n",
1079  nt.n_pntproc,
1080  sz_pntproc,
1081  nt.n_pntproc * sz_pntproc);
1082  printf("n_netcon=%d sz=%ld nbyte=%ld\n", nt.n_netcon, sz_netcon, nt.n_netcon * sz_netcon);
1083  printf("n_weight = %d\n", nt.n_weight);
1084 
1085  printf("%d thread %d total bytes %ld\n", nrnmpi_myid, i, nb_nt);
1086 #endif
1087 
1088  if (detailed_report) {
1089  size_data[0] += nt.ncell;
1090  size_data[1] += nt.end;
1091  size_data[2] += nmech;
1092  size_data[3] += nt._ndata;
1093  size_data[4] += nt._nidata;
1094  size_data[5] += nt._nvdata;
1095  size_data[6] += nt.n_presyn;
1096  size_data[7] += nt.n_input_presyn;
1097  size_data[8] += nt.n_pntproc;
1098  size_data[9] += nt.n_netcon;
1099  size_data[10] += nt.n_weight;
1100  size_data[11] += nb_nt;
1101  }
1102  }
1103 
1104  nbyte += nccnt * sizeof(NetCon*);
1105  nbyte += output_presyn_size();
1106  nbyte += input_presyn_size();
1107 
1108  nbyte += nrnran123_instance_count() * nrnran123_state_size();
1109 
1110 #ifdef DEBUG
1111  printf("%d netcon pointers %ld nbyte=%ld\n", nrnmpi_myid, nccnt, nccnt * sizeof(NetCon*));
1112  printf("nrnran123 size=%ld cnt=%ld nbyte=%ld\n",
1113  nrnran123_state_size(),
1115  nrnran123_instance_count() * nrnran123_state_size());
1116  printf("%d total bytes %ld\n", nrnmpi_myid, nbyte);
1117 #endif
1118  if (detailed_report) {
1119  size_data[12] = nbyte;
1120 #if NRNMPI
1121  if (corenrn_param.mpi_enable) {
1122  // last arg is op type where 1 is sum, 2 is max and any other value is min
1123  nrnmpi_long_allreduce_vec(&size_data[0], &global_size_data_sum[0], 13, 1);
1124  nrnmpi_long_allreduce_vec(&size_data[0], &global_size_data_max[0], 13, 2);
1125  nrnmpi_long_allreduce_vec(&size_data[0], &global_size_data_min[0], 13, 3);
1126  for (int i = 0; i < 13; i++) {
1127  global_size_data_avg[i] = global_size_data_sum[i] / float(nrnmpi_numprocs);
1128  }
1129  } else
1130 #endif
1131  {
1132  global_size_data_max = size_data;
1133  global_size_data_min = size_data;
1134  global_size_data_avg.assign(size_data.cbegin(), size_data.cend());
1135  }
1136  // now print the collected data:
1137  if (nrnmpi_myid == 0) {
1138  printf("Memory size information for all NrnThreads per rank\n");
1139  printf("------------------------------------------------------------------\n");
1140  printf("%22s %12s %12s %12s\n", "field", "min", "max", "avg");
1141  printf("%22s %12ld %12ld %15.2f\n",
1142  "n_cell",
1143  global_size_data_min[0],
1144  global_size_data_max[0],
1145  global_size_data_avg[0]);
1146  printf("%22s %12ld %12ld %15.2f\n",
1147  "n_compartment",
1148  global_size_data_min[1],
1149  global_size_data_max[1],
1150  global_size_data_avg[1]);
1151  printf("%22s %12ld %12ld %15.2f\n",
1152  "n_mechanism",
1153  global_size_data_min[2],
1154  global_size_data_max[2],
1155  global_size_data_avg[2]);
1156  printf("%22s %12ld %12ld %15.2f\n",
1157  "_ndata",
1158  global_size_data_min[3],
1159  global_size_data_max[3],
1160  global_size_data_avg[3]);
1161  printf("%22s %12ld %12ld %15.2f\n",
1162  "_nidata",
1163  global_size_data_min[4],
1164  global_size_data_max[4],
1165  global_size_data_avg[4]);
1166  printf("%22s %12ld %12ld %15.2f\n",
1167  "_nvdata",
1168  global_size_data_min[5],
1169  global_size_data_max[5],
1170  global_size_data_avg[5]);
1171  printf("%22s %12ld %12ld %15.2f\n",
1172  "n_presyn",
1173  global_size_data_min[6],
1174  global_size_data_max[6],
1175  global_size_data_avg[6]);
1176  printf("%22s %12ld %12ld %15.2f\n",
1177  "n_presyn (bytes)",
1178  global_size_data_min[6] * sz_presyn,
1179  global_size_data_max[6] * sz_presyn,
1180  global_size_data_avg[6] * sz_presyn);
1181  printf("%22s %12ld %12ld %15.2f\n",
1182  "n_input_presyn",
1183  global_size_data_min[7],
1184  global_size_data_max[7],
1185  global_size_data_avg[7]);
1186  printf("%22s %12ld %12ld %15.2f\n",
1187  "n_input_presyn (bytes)",
1188  global_size_data_min[7] * sz_input_presyn,
1189  global_size_data_max[7] * sz_input_presyn,
1190  global_size_data_avg[7] * sz_input_presyn);
1191  printf("%22s %12ld %12ld %15.2f\n",
1192  "n_pntproc",
1193  global_size_data_min[8],
1194  global_size_data_max[8],
1195  global_size_data_avg[8]);
1196  printf("%22s %12ld %12ld %15.2f\n",
1197  "n_pntproc (bytes)",
1198  global_size_data_min[8] * sz_pntproc,
1199  global_size_data_max[8] * sz_pntproc,
1200  global_size_data_avg[8] * sz_pntproc);
1201  printf("%22s %12ld %12ld %15.2f\n",
1202  "n_netcon",
1203  global_size_data_min[9],
1204  global_size_data_max[9],
1205  global_size_data_avg[9]);
1206  printf("%22s %12ld %12ld %15.2f\n",
1207  "n_netcon (bytes)",
1208  global_size_data_min[9] * sz_netcon,
1209  global_size_data_max[9] * sz_netcon,
1210  global_size_data_avg[9] * sz_netcon);
1211  printf("%22s %12ld %12ld %15.2f\n",
1212  "n_weight",
1213  global_size_data_min[10],
1214  global_size_data_max[10],
1215  global_size_data_avg[10]);
1216  printf("%22s %12ld %12ld %15.2f\n",
1217  "NrnThread (bytes)",
1218  global_size_data_min[11],
1219  global_size_data_max[11],
1220  global_size_data_avg[11]);
1221  printf("%22s %12ld %12ld %15.2f\n",
1222  "model size (bytes)",
1223  global_size_data_min[12],
1224  global_size_data_max[12],
1225  global_size_data_avg[12]);
1226  }
1227  }
1228 
1229 #if NRNMPI
1230  if (corenrn_param.mpi_enable) {
1231  long global_nbyte = 0;
1232  nrnmpi_long_allreduce_vec(&nbyte, &global_nbyte, 1, 1);
1233  nbyte = global_nbyte;
1234  }
1235 #endif
1236 
1237  return nbyte;
1238 }
1239 
1240 } // namespace coreneuron
int * nrn_has_net_event_
Definition: init.cpp:161
auto & get_memb_func(size_t idx)
Definition: coreneuron.hpp:135
void read_direct(int thread_id, const NrnThread &nt)
Definition: phase2.cpp:279
void read_file(FileHandler &F, const NrnThread &nt)
Definition: phase2.cpp:126
void populate(NrnThread &nt, const UserParams &userParams)
Definition: phase2.cpp:963
void read_file(FileHandler &F, NrnThreadMappingInfo *ntmapping)
Definition: phase3.cpp:29
void read_direct(NrnThreadMappingInfo *ntmapping)
Definition: phase3.cpp:48
static Frame * fp
Definition: code.cpp:96
#define cnt
Definition: tqueue.hpp:44
#define weights
Definition: md1redef.h:42
#define i
Definition: md1redef.h:19
#define assert(ex)
Definition: hocassrt.h:24
struct NrnThread NrnThread
Definition: membfunc.h:15
void free_memory(void *pointer)
Definition: memory.h:213
printf
Definition: extdef.h:5
struct NrnThreadMembList NrnThreadMembList
void gap_data_indices_setup(NrnThread *nt)
For now, until conceptualization of the ordering is clear, just replace src setup_info_ indices value...
TransferThreadData * transfer_thread_data_
Definition: partrans.cpp:25
SetupTransferInfo * setup_info_
void gap_mpi_setup(int ngroup)
THIS FILE IS AUTO GENERATED DONT MODIFY IT.
NrnThread * nrn_threads
Definition: multicore.cpp:56
void nrn_threads_create(int n)
Definition: multicore.cpp:102
int soaos2cnrn_index(const std::array< int, 3 > &soaos_indices, const std::vector< int > &array_dims, int padded_node_count, int *permute)
Compute the CoreNEURON index given an SoAoS index.
void nrn_inverse_i_layout(int i, int &icnt, int cnt, int &isz, int sz, int layout)
Definition: nrn_setup.cpp:680
void clear_event_queue()
Definition: cvodestb.cpp:47
void(*)(NrnThread *, Memb_list *, int) mod_f_t
Definition: membfunc.hpp:24
static OMP_Mutex mut
Definition: nrn_setup.cpp:154
size_t output_presyn_size(void)
Approximate count of number of bytes for the gid2out map.
Definition: nrn_setup.cpp:1008
void nrn_cleanup()
Definition: nrn_setup.cpp:725
void write_mech_report()
display global mechanism count
Definition: mech_report.cpp:19
void read_phase1(NrnThread &nt, UserParams &userParams)
Definition: nrn_setup.cpp:923
void * ecalloc_align(size_t n, size_t size, size_t alignment)
size_t input_presyn_size(void)
Definition: nrn_setup.cpp:1020
int nrn_nthread
Definition: multicore.cpp:55
double ** nrn_ion_global_map
void nrn_setup(const char *filesdat, bool is_mapping_needed, CheckPoints &checkPoints, bool run_setup_cleanup, const char *datpath, const char *restore_path, double *mindelay)
Definition: nrn_setup.cpp:409
std::map< int, InputPreSyn * > gid2in
Definition: nrn_setup.cpp:160
bool nrn_have_gaps
variables defined in coreneuron library
Definition: partrans.cpp:21
void nrn_fast_imem_alloc()
Definition: fast_imem.cpp:32
bool use_multisend_
Definition: multisend.cpp:53
size_t memb_list_size(NrnThreadMembList *tml, bool include_data)
Definition: nrn_setup.cpp:974
static void nrn_fatal_error(const char *msg)
Definition: nrnmpi.cpp:31
void nrn_mk_table_check()
Definition: multicore.cpp:133
void netpar_tid_gid2ps(int tid, int gid, PreSyn **ps, InputPreSyn **psi)
Definition: nrn_setup.cpp:233
int Datum
Definition: nrnconf.h:23
std::vector< std::map< int, PreSyn * > > neg_gid2out
Vector of maps for negative presyns.
Definition: nrn_setup.cpp:157
bool use_solve_interleave
Definition: solve_core.cpp:13
std::vector< std::vector< int > > nrnthreads_netcon_negsrcgid_tid
If a nrnthreads_netcon_srcgid is negative, need to determine the thread when in order to use the corr...
Definition: nrn_setup.cpp:170
void read_phase3(NrnThread &nt, UserParams &userParams)
read mapping information for neurons
Definition: nrn_setup.cpp:941
void check_bbcore_write_version(const char *)
Definition: nrnoc_aux.cpp:128
void nrn_read_filesdat(int &ngrp, int *&grp, const char *filesdat)
Definition: nrn_setup.cpp:173
std::vector< int * > nrnthreads_netcon_srcgid
Only for setup vector of netcon source gids.
Definition: nrn_setup.cpp:166
void delete_trajectory_requests(NrnThread &nt)
Definition: nrn_setup.cpp:905
std::array< int, 3 > legacy2soaos_index(int legacy_index, const std::vector< int > &array_dims)
Split a legacy index into the three SoAoS indices.
void nrn_p_construct()
Definition: netcvode.cpp:175
void read_phasegap(NrnThread &nt, UserParams &userParams)
Definition: nrn_setup.cpp:595
CoreNeuron corenrn
Definition: multicore.cpp:53
int nrn_get_mechtype(const char *name)
Get mechanism type by the mechanism name.
Definition: mk_mech.cpp:145
void setup_ThreadData(NrnThread &nt)
Definition: nrn_setup.cpp:577
void nrn_threads_free()
Definition: multicore.cpp:125
double set_mindelay(double maxdelay)
Definition: netpar.cpp:643
double nrn_wtime()
Definition: utils.cpp:22
NrnThreadChkpnt * nrnthread_chkpnt
void nrn_multisend_setup()
std::size_t nrnran123_instance_count()
Definition: nrnran123.cpp:107
void nrn_setup_cleanup()
Clean up.
Definition: nrn_setup.cpp:399
int nrn_soa_padded_size(int cnt, int layout)
calculate size after padding for specific memory layout
double * legacy_index2pointer(int mtype, int index, NrnThread &nt)
Definition: nrn_setup.cpp:643
void nrn_multithread_job(F &&job, Args &&... args)
Definition: multicore.hpp:161
void determine_inputpresyn()
Definition: nrn_setup.cpp:256
void fast_imem_free()
Definition: fast_imem.cpp:21
std::vector< NetCon * > netcon_in_presyn_order_
InputPreSyn.nc_index_ to + InputPreSyn.nc_cnt_ give the NetCon*.
Definition: nrn_setup.cpp:163
size_t model_size(bool detailed_report)
Definition: nrn_setup.cpp:1032
void nrn_multisend_cleanup()
void delete_fornetcon_info(NrnThread &nt)
Definition: nrn_setup.cpp:716
int nrn_ion_global_map_size
std::map< int, PreSyn * > gid2out
Maps for ouput and input presyns.
Definition: nrn_setup.cpp:159
corenrn_parameters corenrn_param
Printing method.
void nrn_cleanup_ion_map()
Cleanup global ion map created during mechanism registration.
Definition: nrn_setup.cpp:707
void read_phase2(NrnThread &nt, UserParams &userParams)
Definition: nrn_setup.cpp:930
void create_interleave_info()
Definition: cellorder.cpp:110
void destroy_interleave_info()
Definition: cellorder.cpp:115
#define nrn_assert(x)
assert()-like macro, independent of NDEBUG status
Definition: nrn_assert.h:33
int corenrn_embedded_nthread
Definition: nrn_setup.cpp:49
void(* nrn2core_all_weights_return_)(std::vector< double * > &weights)
Definition: nrn_setup.cpp:75
void(* nrn2core_trajectory_values_)(int tid, int n_pr, void **vpr, double t)
Definition: nrn_setup.cpp:69
bool corenrn_file_mode
Definition: nrn_setup.cpp:48
void(* nrn2core_get_trajectory_requests_)(int tid, int &bsize, int &n_pr, void **&vpr, int &n_trajec, int *&types, int *&indices, double **&pvars, double **&varrays)
Definition: nrn_setup.cpp:59
bool corenrn_embedded
--> Coreneuron
Definition: nrn_setup.cpp:47
int(* nrn2core_all_spike_vectors_return_)(std::vector< double > &spikevec, std::vector< int > &gidvec)
Definition: nrn_setup.cpp:73
void(* nrn2core_group_ids_)(int *)
Definition: nrn_setup.cpp:51
void(* nrn2core_trajectory_return_)(int tid, int n_pr, int bsize, int vecsz, void **vpr, double t)
Definition: nrn_setup.cpp:71
int const size_t const size_t n
Definition: nrngsl.h:10
s
Definition: multisend.cpp:521
short index
Definition: cabvars.h:11
std::vector< Memb_func > memb_func
Definition: init.cpp:145
#define lock
int sgid_t
Definition: partrans.hpp:20
std::vector< sgid_t > src_sid
Definition: partrans.hpp:95
void(* thread_mem_init_)(ThreadDatum *)
Definition: membfunc.hpp:48
ThreadDatum * _thread
Definition: mechanism.hpp:141
NetSendBuffer_t * _net_send_buffer
Definition: mechanism.hpp:143
NetReceiveBuffer_t * _net_receive_buffer
Definition: mechanism.hpp:142
PreSynHelper * presyns_helper
Definition: multicore.hpp:84
NrnFastImem * nrn_fast_imem
Definition: multicore.hpp:124
size_t * _fornetcon_weight_perm
Definition: multicore.hpp:152
size_t * _fornetcon_perm_indices
Definition: multicore.hpp:150
Memb_list ** _ml_list
Definition: multicore.hpp:81
std::unique_ptr< SummationReportMapping > summation_report_handler_
Definition: multicore.hpp:144
NrnThreadMembList * tml
Definition: multicore.hpp:80
TrajectoryRequests * trajec_requests
Definition: multicore.hpp:146
Point_process * pntprocs
Definition: multicore.hpp:82
Compartment mapping information for NrnThread.
size_t size() const
number of cells
NrnThreadMembList * next
Definition: multicore.hpp:33
This structure is data needed is several part of nrn_setup, phase1 and phase2.
Definition: user_params.hpp:18
const int ngroup
direct memory mode with neuron, do not open files Number of local cell groups
Definition: user_params.hpp:33
const int *const gidgroups
Array of cell group numbers (indices)
Definition: user_params.hpp:35
std::vector< FileHandler > file_reader
Definition: user_params.hpp:40
bool model_stats
Print version and exit.
unsigned cell_interleave_permute
Spike Compression.
bool mpi_enable
Initialization seed for random number generator (int)
The basic problem is to copy sources to targets.
Definition: partrans.hpp:78
Project version information.
Definition: config.h:26