NEURON
phase2.cpp
Go to the documentation of this file.
1 /*
2 # =============================================================================
3 # Copyright (c) 2016 - 2022 Blue Brain Project/EPFL
4 #
5 # See top-level LICENSE file for details.
6 # =============================================================================
7 */
8 
9 #include "nrnoc/ion_semantics.h"
10 #include "coreneuron/io/phase2.hpp"
22 
23 #if defined(_OPENMP)
24 #include <omp.h>
25 #endif
26 
27 int (*nrn2core_get_dat2_1_)(int tid,
28  int& n_real_cell,
29  int& ngid,
30  int& n_real_gid,
31  int& nnode,
32  int& ndiam,
33  int& nmech,
34  int*& tml_index,
35  int*& ml_nodecount,
36  int& nidata,
37  int& nvdata,
38  int& nweight);
39 
40 int (*nrn2core_get_dat2_2_)(int tid,
41  int*& v_parent_index,
42  double*& a,
43  double*& b,
44  double*& area,
45  double*& v,
46  double*& diamvec);
47 
48 int (*nrn2core_get_dat2_mech_)(int tid,
49  size_t i,
50  int dsz_inst,
51  int*& nodeindices,
52  double*& data,
53  int*& pdata,
54  std::vector<uint32_t>& nmodlrandom,
55  std::vector<int>& pointer2type);
56 
57 int (*nrn2core_get_dat2_3_)(int tid,
58  int nweight,
59  int*& output_vindex,
60  double*& output_threshold,
61  int*& netcon_pnttype,
62  int*& netcon_pntindex,
63  double*& weights,
64  double*& delays);
65 
66 int (*nrn2core_get_dat2_corepointer_)(int tid, int& n);
67 
69  int type,
70  int& icnt,
71  int& dcnt,
72  int*& iarray,
73  double*& darray);
74 
75 int (*nrn2core_get_dat2_vecplay_)(int tid, std::vector<int>& indices);
76 
78  int i,
79  int& vptype,
80  int& mtype,
81  int& ix,
82  int& sz,
83  double*& yvec,
84  double*& tvec,
85  int& last_index,
86  int& discon_index,
87  int& ubound_index);
88 
89 namespace coreneuron {
90 template <typename T>
91 void mech_data_layout_transform(T* data, int cnt, const std::vector<int>& array_dims, int layout) {
92  if (layout == Layout::AoS) {
93  throw std::runtime_error("AoS memory layout not implemented.");
94  }
95 
96  int n_vars = array_dims.size();
97  int row_width = std::accumulate(array_dims.begin(), array_dims.end(), 0);
98  int padded_cnt = nrn_soa_padded_size(cnt, layout);
99 
100  std::vector<T> tmp(padded_cnt * row_width);
101  std::copy(data, data + cnt * row_width, tmp.begin());
102 
103  size_t offset_var = 0;
104  for (size_t i_var = 0; i_var < n_vars; ++i_var) {
105  size_t K = array_dims[i_var];
106 
107  for (size_t i = 0; i < cnt; ++i) {
108  for (size_t k = 0; k < K; ++k) {
109  size_t i_dst = padded_cnt * offset_var + i * K + k;
110  size_t i_src = i * row_width + offset_var + k;
111  data[i_dst] = tmp[i_src];
112  }
113  }
114 
115  offset_var += K;
116  }
117 }
118 
119 template <typename T>
120 inline void mech_data_layout_transform(T* data, int cnt, int sz, int layout) {
121  std::vector<int> array_dims(sz, 1);
122  mech_data_layout_transform(data, cnt, array_dims, layout);
123 }
124 
125 
127  n_real_cell = F.read_int();
128  n_output = F.read_int();
129  n_real_output = F.read_int();
130  n_node = F.read_int();
131  n_diam = F.read_int();
132  n_mech = F.read_int();
133  mech_types = std::vector<int>(n_mech, 0);
134  nodecounts = std::vector<int>(n_mech, 0);
135  for (int i = 0; i < n_mech; ++i) {
136  mech_types[i] = F.read_int();
137  nodecounts[i] = F.read_int();
138  }
139 
140  // check mechanism compatibility before reading data
141  check_mechanism();
142 
143  n_idata = F.read_int();
144  n_vdata = F.read_int();
145  int n_weight = F.read_int();
146  v_parent_index = (int*) ecalloc_align(n_node, sizeof(int));
148 
149  int n_data_padded = nrn_soa_padded_size(n_node, SOA_LAYOUT);
150  {
151  { // Compute size of _data and allocate
152  int n_data = 6 * n_data_padded;
153  if (n_diam > 0) {
154  n_data += n_data_padded;
155  }
156  for (int i = 0; i < n_mech; ++i) {
157  int layout = corenrn.get_mech_data_layout()[mech_types[i]];
158  int n = nodecounts[i];
160  n_data = nrn_soa_byte_align(n_data);
161  n_data += nrn_soa_padded_size(n, layout) * sz;
162  }
163  _data = (double*) ecalloc_align(n_data, sizeof(double));
164  }
165  F.read_array<double>(_data + 2 * n_data_padded, n_node);
166  F.read_array<double>(_data + 3 * n_data_padded, n_node);
167  F.read_array<double>(_data + 5 * n_data_padded, n_node);
168  F.read_array<double>(_data + 4 * n_data_padded, n_node);
169  if (n_diam > 0) {
170  F.read_array<double>(_data + 6 * n_data_padded, n_node);
171  }
172  }
173 
174  size_t offset = 6 * n_data_padded;
175  if (n_diam > 0) {
176  offset += n_data_padded;
177  }
178  for (int i = 0; i < n_mech; ++i) {
179  int layout = corenrn.get_mech_data_layout()[mech_types[i]];
180  int n = nodecounts[i];
183  offset = nrn_soa_byte_align(offset);
184  std::vector<int> nodeindices;
186  nodeindices = F.read_vector<int>(n);
187  }
188  F.read_array<double>(_data + offset, sz * n);
189  offset += nrn_soa_padded_size(n, layout) * sz;
190  std::vector<int> pdata;
191  if (dsz > 0) {
192  pdata = F.read_vector<int>(dsz * n);
193  }
194  tmls.emplace_back(TML{nodeindices, pdata, mech_types[i], {}, {}});
195  if (dsz > 0) {
196  int sz = F.read_int();
197  if (sz) {
198  auto& p2t = tmls.back().pointer2type;
199  p2t = F.read_vector<int>(sz);
200  }
201 
202  // nmodlrandom
203  sz = F.read_int();
204  if (sz) {
205  auto& nmodlrandom = tmls.back().nmodlrandom;
206  nmodlrandom = F.read_vector<uint32_t>(sz);
207  }
208  }
209  }
210  output_vindex = F.read_vector<int>(nt.n_presyn);
212  pnttype = F.read_vector<int>(nt.n_netcon);
213  pntindex = F.read_vector<int>(nt.n_netcon);
214  weights = F.read_vector<double>(n_weight);
215  delay = F.read_vector<double>(nt.n_netcon);
217 
218  for (int i = 0; i < n_mech; ++i) {
219  if (!corenrn.get_bbcore_read()[mech_types[i]]) {
220  continue;
221  }
222  tmls[i].type = F.read_int();
223  int icnt = F.read_int();
224  int dcnt = F.read_int();
225  if (icnt > 0) {
226  tmls[i].iArray = F.read_vector<int>(icnt);
227  }
228  if (dcnt > 0) {
229  tmls[i].dArray = F.read_vector<double>(dcnt);
230  }
231  }
232 
233  int n_vec_play_continuous = F.read_int();
234  vec_play_continuous.reserve(n_vec_play_continuous);
235  for (int i = 0; i < n_vec_play_continuous; ++i) {
236  VecPlayContinuous_ item;
237  item.vtype = F.read_int();
238  item.mtype = F.read_int();
239  item.ix = F.read_int();
240  int sz = F.read_int();
241  item.yvec = IvocVect(sz);
242  item.tvec = IvocVect(sz);
243  F.read_array<double>(item.yvec.data(), sz);
244  F.read_array<double>(item.tvec.data(), sz);
245  vec_play_continuous.push_back(std::move(item));
246  }
247 
248  // store current checkpoint state to continue reading mapping
249  // The checkpoint numbering in phase 3 is a continuing of phase 2, and so will be restored
250  F.record_checkpoint();
251 
252  if (F.eof())
253  return;
254 
255  nrn_assert(F.read_int() == n_vec_play_continuous);
256 
257  for (int i = 0; i < n_vec_play_continuous; ++i) {
258  auto& vecPlay = vec_play_continuous[i];
259  vecPlay.last_index = F.read_int();
260  vecPlay.discon_index = F.read_int();
261  vecPlay.ubound_index = F.read_int();
262  }
263 
264  patstim_index = F.read_int();
265 
266  nrn_assert(F.read_int() == -1);
267 
268  for (int i = 0; i < nt.n_presyn; ++i) {
269  preSynConditionEventFlags.push_back(F.read_int());
270  }
271 
272  nrn_assert(F.read_int() == -1);
273  restore_events(F);
274 
275  nrn_assert(F.read_int() == -1);
276  restore_events(F);
277 }
278 
279 void Phase2::read_direct(int thread_id, const NrnThread& nt) {
280  int* types_ = nullptr;
281  int* nodecounts_ = nullptr;
282  int n_weight;
283  (*nrn2core_get_dat2_1_)(thread_id,
284  n_real_cell,
285  n_output,
287  n_node,
288  n_diam,
289  n_mech,
290  types_,
291  nodecounts_,
292  n_idata,
293  n_vdata,
294  n_weight);
295  mech_types = std::vector<int>(types_, types_ + n_mech);
296  delete[] types_;
297 
298  nodecounts = std::vector<int>(nodecounts_, nodecounts_ + n_mech);
299  delete[] nodecounts_;
300 
301  check_mechanism();
302 
303  // TODO: fix it in the future
304  int n_data_padded = nrn_soa_padded_size(n_node, SOA_LAYOUT);
305  int n_data = 6 * n_data_padded;
306  if (n_diam > 0) {
307  n_data += n_data_padded;
308  }
309  for (int i = 0; i < n_mech; ++i) {
310  int layout = corenrn.get_mech_data_layout()[mech_types[i]];
311  int n = nodecounts[i];
313  n_data = nrn_soa_byte_align(n_data);
314  n_data += nrn_soa_padded_size(n, layout) * sz;
315  }
316  _data = (double*) ecalloc_align(n_data, sizeof(double));
317 
318  v_parent_index = (int*) ecalloc_align(n_node, sizeof(int));
319  double* actual_a = _data + 2 * n_data_padded;
320  double* actual_b = _data + 3 * n_data_padded;
321  double* actual_v = _data + 4 * n_data_padded;
322  double* actual_area = _data + 5 * n_data_padded;
323  double* actual_diam = n_diam > 0 ? _data + 6 * n_data_padded : nullptr;
324  (*nrn2core_get_dat2_2_)(
325  thread_id, v_parent_index, actual_a, actual_b, actual_area, actual_v, actual_diam);
326 
327  tmls.resize(n_mech);
328 
329  auto& param_sizes = corenrn.get_prop_param_size();
330  auto& dparam_sizes = corenrn.get_prop_dparam_size();
331  int dsz_inst = 0;
332  size_t offset = 6 * n_data_padded;
333  if (n_diam > 0) {
334  offset += n_data_padded;
335  }
336  for (int i = 0; i < n_mech; ++i) {
337  auto& tml = tmls[i];
338  int type = mech_types[i];
339  int layout = corenrn.get_mech_data_layout()[type];
340  offset = nrn_soa_byte_align(offset);
341 
342  tml.type = type;
343  // artificial cell don't use nodeindices
344  if (!corenrn.get_is_artificial()[type]) {
345  tml.nodeindices.resize(nodecounts[i]);
346  }
347  tml.pdata.resize(nodecounts[i] * dparam_sizes[type]);
348 
349  int* nodeindices_ = nullptr;
350  double* data_ = _data + offset;
351  int* pdata_ = const_cast<int*>(tml.pdata.data());
352 
353  // nmodlrandom, receives:
354  // number of random variables
355  // dparam index (neuron side) of each random variable
356  // 5 uint32 for each var of each instance
357  // id1, id2, id3, seq, uint32_t(which)
358  // all instances of ranvar1 first, then all instances of ranvar2, etc.
359  auto& nmodlrandom = tml.nmodlrandom;
360 
361  (*nrn2core_get_dat2_mech_)(thread_id,
362  i,
363  dparam_sizes[type] > 0 ? dsz_inst : 0,
364  nodeindices_,
365  data_,
366  pdata_,
367  nmodlrandom,
368  tml.pointer2type);
369 
370  if (dparam_sizes[type] > 0) {
371  dsz_inst++;
372  }
373  offset += nrn_soa_padded_size(nodecounts[i], layout) * param_sizes[type];
374  if (nodeindices_) {
375  std::copy(nodeindices_, nodeindices_ + nodecounts[i], tml.nodeindices.data());
376  free(nodeindices_); // not free_memory because this is allocated by NEURON?
377  }
378  if (corenrn.get_is_artificial()[type]) {
379  assert(nodeindices_ == nullptr);
380  }
381  }
382 
383  int* output_vindex_ = nullptr;
384  double* output_threshold_ = nullptr;
385  int* pnttype_ = nullptr;
386  int* pntindex_ = nullptr;
387  double* weight_ = nullptr;
388  double* delay_ = nullptr;
389  (*nrn2core_get_dat2_3_)(thread_id,
390  n_weight,
391  output_vindex_,
392  output_threshold_,
393  pnttype_,
394  pntindex_,
395  weight_,
396  delay_);
397 
398  output_vindex = std::vector<int>(output_vindex_, output_vindex_ + nt.n_presyn);
399  delete[] output_vindex_;
400 
401  output_threshold = std::vector<double>(output_threshold_, output_threshold_ + n_real_output);
402  delete[] output_threshold_;
403 
404  int n_netcon = nt.n_netcon;
405  pnttype = std::vector<int>(pnttype_, pnttype_ + n_netcon);
406  delete[] pnttype_;
407 
408  pntindex = std::vector<int>(pntindex_, pntindex_ + n_netcon);
409  delete[] pntindex_;
410 
411  weights = std::vector<double>(weight_, weight_ + n_weight);
412  delete[] weight_;
413 
414  delay = std::vector<double>(delay_, delay_ + n_netcon);
415  delete[] delay_;
416 
417  (*nrn2core_get_dat2_corepointer_)(nt.id, num_point_process);
418 
419  for (int i = 0; i < n_mech; ++i) {
420  // not all mod files have BBCOREPOINTER data to read
421  if (!corenrn.get_bbcore_read()[mech_types[i]]) {
422  continue;
423  }
424  int icnt;
425  int* iArray_ = nullptr;
426  int dcnt;
427  double* dArray_ = nullptr;
428  (*nrn2core_get_dat2_corepointer_mech_)(nt.id, tmls[i].type, icnt, dcnt, iArray_, dArray_);
429  tmls[i].iArray.resize(icnt);
430  std::copy(iArray_, iArray_ + icnt, tmls[i].iArray.begin());
431  delete[] iArray_;
432 
433  tmls[i].dArray.resize(dcnt);
434  std::copy(dArray_, dArray_ + dcnt, tmls[i].dArray.begin());
435  delete[] dArray_;
436  }
437 
438  // Get from NEURON, the VecPlayContinuous indices in
439  // NetCvode::fixed_play_ for this thread.
440  std::vector<int> indices_vec_play_continuous;
441  (*nrn2core_get_dat2_vecplay_)(thread_id, indices_vec_play_continuous);
442 
443  // i is an index into NEURON's NetCvode::fixed_play_ for this thread.
444  for (auto i: indices_vec_play_continuous) {
445  VecPlayContinuous_ item;
446  // yvec_ and tvec_ are not deleted as that space is within
447  // NEURON Vector
448  double *yvec_, *tvec_;
449  int sz;
450  (*nrn2core_get_dat2_vecplay_inst_)(thread_id,
451  i,
452  item.vtype,
453  item.mtype,
454  item.ix,
455  sz,
456  yvec_,
457  tvec_,
458  item.last_index,
459  item.discon_index,
460  item.ubound_index);
461  item.yvec = IvocVect(sz);
462  item.tvec = IvocVect(sz);
463  std::copy(yvec_, yvec_ + sz, item.yvec.data());
464  std::copy(tvec_, tvec_ + sz, item.tvec.data());
465  vec_play_continuous.push_back(std::move(item));
466  }
467 }
468 
469 /// Check if MOD file used between NEURON and CoreNEURON is same
471  int diff_mech_count = 0;
472  for (int i = 0; i < n_mech; ++i) {
473  if (std::any_of(corenrn.get_different_mechanism_type().begin(),
475  [&](int e) { return e == mech_types[i]; })) {
476  if (nrnmpi_myid == 0) {
477  printf("Error: %s is a different MOD file than used by NEURON!\n",
479  }
480  diff_mech_count++;
481  }
482  }
483 
484  if (diff_mech_count > 0) {
485  if (nrnmpi_myid == 0) {
486  printf(
487  "Error : NEURON and CoreNEURON must use same mod files for compatibility, %d "
488  "different mod file(s) found. Re-compile special and special-core!\n",
489  diff_mech_count);
490  nrn_abort(1);
491  }
492  }
493 }
494 
495 /// Perform in memory transformation between AoS<>SoA for integer data
497  int nodecount,
498  int* pdata,
499  int i,
500  int dparam_size,
501  int layout,
502  int n_node_) {
503  for (int iml = 0; iml < nodecount; ++iml) {
504  int* pd = pdata + nrn_i_layout(iml, nodecount, i, dparam_size, layout);
505  int ix = *pd; // relative to beginning of _actual_*
506  nrn_assert((ix >= 0) && (ix < n_node_));
507  *pd = elem0 + ix; // relative to nt._data
508  }
509 }
510 
511 void Phase2::set_net_send_buffer(Memb_list** ml_list, const std::vector<int>& pnt_offset) {
512  // NetReceiveBuffering
513  for (auto& net_buf_receive: corenrn.get_net_buf_receive()) {
514  int type = net_buf_receive.second;
515  // Does this thread have this type.
516  Memb_list* ml = ml_list[type];
517  if (ml) { // needs a NetReceiveBuffer
518  NetReceiveBuffer_t* nrb =
521  ml->_net_receive_buffer = nrb;
522  nrb->_pnt_offset = pnt_offset[type];
523 
524  // begin with a size equal to the number of instances, or at least 8
525  nrb->_size = std::max(8, ml->nodecount);
526  nrb->_pnt_index = (int*) ecalloc_align(nrb->_size, sizeof(int));
527  nrb->_displ = (int*) ecalloc_align(nrb->_size + 1, sizeof(int));
528  nrb->_nrb_index = (int*) ecalloc_align(nrb->_size, sizeof(int));
529  nrb->_weight_index = (int*) ecalloc_align(nrb->_size, sizeof(int));
530  nrb->_nrb_t = (double*) ecalloc_align(nrb->_size, sizeof(double));
531  nrb->_nrb_flag = (double*) ecalloc_align(nrb->_size, sizeof(double));
532  }
533  }
534 
535  // NetSendBuffering
536  for (int type: corenrn.get_net_buf_send_type()) {
537  // Does this thread have this type.
538  Memb_list* ml = ml_list[type];
539  if (ml) { // needs a NetSendBuffer
540  assert(!ml->_net_send_buffer);
541  // begin with a size equal to twice number of instances
542  NetSendBuffer_t* nsb = new NetSendBuffer_t(ml->nodecount * 2);
543  ml->_net_send_buffer = nsb;
544  }
545  }
546 }
547 
549  int type;
550  while ((type = F.read_int()) != 0) {
551  double time;
552  F.read_array(&time, 1);
553  switch (type) {
554  case NetConType: {
555  auto event = std::make_shared<NetConType_>();
556  event->time = time;
557  event->netcon_index = F.read_int();
558  events.emplace_back(type, event);
559  break;
560  }
561  case SelfEventType: {
562  auto event = std::make_shared<SelfEventType_>();
563  event->time = time;
564  event->target_type = F.read_int();
565  event->point_proc_instance = F.read_int();
566  event->target_instance = F.read_int();
567  F.read_array(&event->flag, 1);
568  event->movable = F.read_int();
569  event->weight_index = F.read_int();
570  events.emplace_back(type, event);
571  break;
572  }
573  case PreSynType: {
574  auto event = std::make_shared<PreSynType_>();
575  event->time = time;
576  event->presyn_index = F.read_int();
577  events.emplace_back(type, event);
578  break;
579  }
580  case NetParEventType: {
581  auto event = std::make_shared<NetParEvent_>();
582  event->time = time;
583  events.emplace_back(type, event);
584  break;
585  }
586  case PlayRecordEventType: {
587  auto event = std::make_shared<PlayRecordEventType_>();
588  event->time = time;
589  event->play_record_type = F.read_int();
590  if (event->play_record_type == VecPlayContinuousType) {
591  event->vecplay_index = F.read_int();
592  events.emplace_back(type, event);
593  } else {
594  nrn_assert(0);
595  }
596  break;
597  }
598  default: {
599  nrn_assert(0);
600  break;
601  }
602  }
603  }
604 }
605 
606 void Phase2::fill_before_after_lists(NrnThread& nt, const std::vector<Memb_func>& memb_func) {
607  /// Fill the BA lists
608  std::vector<BAMech*> before_after_map(memb_func.size());
609  for (int i = 0; i < BEFORE_AFTER_SIZE; ++i) {
610  for (size_t ii = 0; ii < memb_func.size(); ++ii) {
611  before_after_map[ii] = nullptr;
612  }
613  // Save first before-after block only. In case of multiple before-after blocks with the
614  // same mech type, we will get subsequent ones using linked list below.
615  for (auto bam = corenrn.get_bamech()[i]; bam; bam = bam->next) {
616  if (!before_after_map[bam->type]) {
617  before_after_map[bam->type] = bam;
618  }
619  }
620  // necessary to keep in order wrt multiple BAMech with same mech type
621  NrnThreadBAList** ptbl = nt.tbl + i;
622  for (auto tml = nt.tml; tml; tml = tml->next) {
623  if (before_after_map[tml->index]) {
624  int mtype = tml->index;
625  for (auto bam = before_after_map[mtype]; bam && bam->type == mtype;
626  bam = bam->next) {
627  auto tbl = (NrnThreadBAList*) emalloc(sizeof(NrnThreadBAList));
628  *ptbl = tbl;
629  tbl->next = nullptr;
630  tbl->bam = bam;
631  tbl->ml = tml->ml;
632  ptbl = &(tbl->next);
633  }
634  }
635  }
636  }
637 }
638 
639 void Phase2::pdata_relocation(const NrnThread& nt, const std::vector<Memb_func>& memb_func) {
640  // Some pdata may index into data which has been reordered from AoS to
641  // SoA. The four possibilities are if semantics is -1 (area), -5 (pointer),
642  // -9 (diam), // or 0-999 (ion variables).
643  // Note that pdata has a layout and the // type block in nt.data into which
644  // it indexes, has a layout.
645 
646  // For faster search of tmls[i].type == type, use a map.
647  // (perhaps would be better to replace tmls so that we can use tmls[type].
648  std::map<int, size_t> type2itml;
649  for (size_t i = 0; i < tmls.size(); ++i) {
650  if (tmls[i].pointer2type.size()) {
651  type2itml[tmls[i].type] = i;
652  }
653  }
654 
655  for (auto tml = nt.tml; tml; tml = tml->next) {
656  int type = tml->index;
657  int layout = corenrn.get_mech_data_layout()[type];
658  int* pdata = tml->ml->pdata;
659  int cnt = tml->ml->nodecount;
660  int szdp = corenrn.get_prop_dparam_size()[type];
661  int* semantics = memb_func[type].dparam_semantics;
662 
663  // compute only for ARTIFICIAL_CELL (has useful area pointer with semantics=-1)
664  if (!corenrn.get_is_artificial()[type]) {
665  if (szdp) {
666  if (!semantics)
667  continue; // temporary for HDFReport, Binreport which will be skipped in
668  // bbcore_write of HBPNeuron
669  nrn_assert(semantics);
670  }
671 
672  for (int i = 0; i < szdp; ++i) {
673  int s = semantics[i];
674  switch (s) {
675  case -1: // area
677  nt._actual_area - nt._data, cnt, pdata, i, szdp, layout, nt.end);
678  break;
679  case -9: // diam
681  nt._actual_diam - nt._data, cnt, pdata, i, szdp, layout, nt.end);
682  break;
683  case -5: // pointer assumes a pointer to membrane voltage
684  // or mechanism data in this thread. The value of the
685  // pointer on the NEURON side was analyzed by
686  // nrn_dblpntr2nrncore which returned the
687  // mechanism index and type. At this moment the index
688  // is in pdata and the type is in tmls[type].pointer2type.
689  // However the latter order is according to the nested
690  // iteration for nodecount { for szdp {}}
691  // Also the nodecount POINTER instances of mechanism
692  // might possibly point to differnt range variables.
693  // Therefore it is not possible to use transform_int_data
694  // and the transform must be done one at a time.
695  // So we do nothing here and separately iterate
696  // after this loop instead of the former voltage only
697  /**
698  transform_int_data(
699  nt._actual_v - nt._data, cnt, pdata, i, szdp, layout, nt.end);
700  **/
701  break;
702  default:
703  if (nrn_semantics_is_ion(s)) { // ion
704  int etype = nrn_semantics_ion_type(s);
705  /* if ion is SoA, must recalculate pdata values */
706  /* if ion is AoS, have to deal with offset */
707  Memb_list* eml = nt._ml_list[etype];
708  int edata0 = eml->data - nt._data;
709  int ecnt = eml->nodecount;
710  int ecnt_padded = nrn_soa_padded_size(ecnt, Layout::SoA);
711  int esz = corenrn.get_prop_param_size()[etype];
712  const std::vector<int>& array_dims = corenrn.get_array_dims()[etype];
713  for (int iml = 0; iml < cnt; ++iml) {
714  int* pd = pdata + nrn_i_layout(iml, cnt, i, szdp, layout);
715  int legacy_index = *pd; // relative to the ion data
716  nrn_assert((legacy_index >= 0) && (legacy_index < ecnt * esz));
717 
718  auto soaos_index = legacy2soaos_index(legacy_index, array_dims);
719  *pd = edata0 +
720  soaos2cnrn_index(soaos_index, array_dims, ecnt_padded, nullptr);
721  }
722  }
723  }
724  }
725  // Handle case -5 POINTER transformation (see comment above)
726  auto search = type2itml.find(type);
727  if (search != type2itml.end()) {
728  auto& ptypes = tmls[type2itml[type]].pointer2type;
729  assert(ptypes.size());
730  size_t iptype = 0;
731  for (int iml = 0; iml < cnt; ++iml) {
732  for (int i = 0; i < szdp; ++i) {
733  if (semantics[i] == -5) { // POINTER
734  int* pd = pdata + nrn_i_layout(iml, cnt, i, szdp, layout);
735  int ix = *pd; // relative to elem0
736  int ptype = ptypes[iptype++];
737  if (ptype == voltage) {
738  nrn_assert((ix >= 0) && (ix < nt.end));
739  int elem0 = nt._actual_v - nt._data;
740  *pd = elem0 + ix;
741  } else {
742  Memb_list* pml = nt._ml_list[ptype];
743  int pcnt = pml->nodecount;
744  int pcnt_padded = nrn_soa_padded_size(pcnt, Layout::SoA);
745  int psz = corenrn.get_prop_param_size()[ptype];
746  const std::vector<int>& array_dims =
747  corenrn.get_array_dims()[ptype];
748  nrn_assert((ix >= 0) && (ix < pcnt * psz));
749 
750  int elem0 = pml->data - nt._data;
751  auto soaos_index = legacy2soaos_index(ix, array_dims);
752  *pd =
753  elem0 +
754  soaos2cnrn_index(soaos_index, array_dims, pcnt_padded, nullptr);
755  }
756  }
757  }
758  }
759  ptypes.clear();
760  }
761  }
762  }
763 }
764 
765 void Phase2::set_dependencies(const NrnThread& nt, const std::vector<Memb_func>& memb_func) {
766  /* here we setup the mechanism dependencies. if there is a mechanism dependency
767  * then we allocate an array for tml->dependencies otherwise set it to nullptr.
768  * In order to find out the "real" dependencies i.e. dependent mechanism
769  * exist at the same compartment, we compare the nodeindices of mechanisms
770  * returned by nrn_mech_depend.
771  */
772 
773  /* temporary array for dependencies */
774  int* mech_deps = (int*) ecalloc(memb_func.size(), sizeof(int));
775 
776  for (auto tml = nt.tml; tml; tml = tml->next) {
777  /* initialize to null */
778  tml->dependencies = nullptr;
779  tml->ndependencies = 0;
780 
781  /* get dependencies from the models */
782  int deps_cnt = nrn_mech_depend(tml->index, mech_deps);
783 
784  /* if dependencies, setup dependency array */
785  if (deps_cnt) {
786  /* store "real" dependencies in the vector */
787  std::vector<int> actual_mech_deps;
788 
789  Memb_list* ml = tml->ml;
790  int* nodeindices = ml->nodeindices;
791 
792  /* iterate over dependencies */
793  for (int j = 0; j < deps_cnt; j++) {
794  /* memb_list of dependency mechanism */
795  Memb_list* dml = nt._ml_list[mech_deps[j]];
796 
797  /* dependency mechanism may not exist in the model */
798  if (!dml)
799  continue;
800 
801  /* take nodeindices for comparison */
802  int* dnodeindices = dml->nodeindices;
803 
804  /* set_intersection function needs temp vector to push the common values */
805  std::vector<int> node_intersection;
806 
807  /* make sure they have non-zero nodes and find their intersection */
808  if ((ml->nodecount > 0) && (dml->nodecount > 0)) {
809  std::set_intersection(nodeindices,
810  nodeindices + ml->nodecount,
811  dnodeindices,
812  dnodeindices + dml->nodecount,
813  std::back_inserter(node_intersection));
814  }
815 
816  /* if they intersect in the nodeindices, it's real dependency */
817  if (!node_intersection.empty()) {
818  actual_mech_deps.push_back(mech_deps[j]);
819  }
820  }
821 
822  /* copy actual_mech_deps to dependencies */
823  if (!actual_mech_deps.empty()) {
824  tml->ndependencies = actual_mech_deps.size();
825  tml->dependencies = (int*) ecalloc(actual_mech_deps.size(), sizeof(int));
826  std::copy(actual_mech_deps.begin(), actual_mech_deps.end(), tml->dependencies);
827  }
828  }
829  }
830 
831  /* free temp dependency array */
832  free(mech_deps);
833 }
834 
835 void Phase2::handle_weights(NrnThread& nt, int n_netcon, NrnThreadChkpnt& ntc) {
836  nt.n_weight = weights.size();
837  // weights in netcons order in groups defined by Point_process target type.
838  nt.weights = (double*) ecalloc_align(nt.n_weight, sizeof(double));
839  std::copy(weights.begin(), weights.end(), nt.weights);
840 
841  int iw = 0;
842  for (int i = 0; i < n_netcon; ++i) {
843  NetCon& nc = nt.netcons[i];
844  nc.u.weight_index_ = iw;
845  if (pnttype[i] != 0) {
847  } else {
848  iw += 1;
849  }
850  }
851  assert(iw == nt.n_weight);
852 
853  // Nontrivial if FOR_NETCON in use by some mechanisms
855 
856 
857 #if CHKPNTDEBUG
858  ntc.delay = new double[n_netcon];
859  memcpy(ntc.delay, delay.data(), n_netcon * sizeof(double));
860 #endif
861  for (int i = 0; i < n_netcon; ++i) {
862  NetCon& nc = nt.netcons[i];
863  nc.delay_ = delay[i];
864  }
865 }
866 
868  const std::vector<Memb_func>& memb_func,
869  NrnThreadChkpnt& ntc) {
870  // BBCOREPOINTER information
871 #if CHKPNTDEBUG
872  ntc.nbcp = num_point_process;
873  ntc.bcpicnt = new int[n_mech];
874  ntc.bcpdcnt = new int[n_mech];
875  ntc.bcptype = new int[n_mech];
876  size_t point_proc_id = 0;
877 #endif
878  for (int i = 0; i < n_mech; ++i) {
879  int type = mech_types[i];
880  if (!corenrn.get_bbcore_read()[type]) {
881  continue;
882  }
883  type = tmls[i].type; // This is not an error, but it has to be fixed I think
884 #if CHKPNTDEBUG
885  ntc.bcptype[point_proc_id] = type;
886  ntc.bcpicnt[point_proc_id] = tmls[i].iArray.size();
887  ntc.bcpdcnt[point_proc_id] = tmls[i].dArray.size();
888  point_proc_id++;
889 #endif
890  int ik = 0;
891  int dk = 0;
892  Memb_list* ml = nt._ml_list[type];
893  int dsz = corenrn.get_prop_param_size()[type];
894  int pdsz = corenrn.get_prop_dparam_size()[type];
895  int cntml = ml->nodecount;
896  int layout = corenrn.get_mech_data_layout()[type];
897  for (int j = 0; j < cntml; ++j) {
898  int jp = j;
899  if (ml->_permute) {
900  jp = ml->_permute[j];
901  }
902  double* d = ml->data;
903  Datum* pd = ml->pdata;
904  d += nrn_i_layout(jp, cntml, 0, dsz, layout);
905  pd += nrn_i_layout(jp, cntml, 0, pdsz, layout);
906  int aln_cntml = nrn_soa_padded_size(cntml, layout);
907  (*corenrn.get_bbcore_read()[type])(tmls[i].dArray.data(),
908  tmls[i].iArray.data(),
909  &dk,
910  &ik,
911  0,
912  aln_cntml,
913  d,
914  pd,
915  ml->_thread,
916  &nt,
917  ml,
918  0.0);
919  }
920  assert(dk == static_cast<int>(tmls[i].dArray.size()));
921  assert(ik == static_cast<int>(tmls[i].iArray.size()));
922  }
923 }
924 
926  // VecPlayContinuous instances
927  // No attempt at memory efficiency
928  nt.n_vecplay = vec_play_continuous.size();
929  if (nt.n_vecplay) {
930  nt._vecplay = new void*[nt.n_vecplay];
931  } else {
932  nt._vecplay = nullptr;
933  }
934 #if CHKPNTDEBUG
935  ntc.vecplay_ix = new int[nt.n_vecplay];
936  ntc.vtype = new int[nt.n_vecplay];
937  ntc.mtype = new int[nt.n_vecplay];
938 #endif
939  for (int i = 0; i < nt.n_vecplay; ++i) {
940  auto& vecPlay = vec_play_continuous[i];
941  nrn_assert(vecPlay.vtype == VecPlayContinuousType);
942 #if CHKPNTDEBUG
943  ntc.vtype[i] = vecPlay.vtype;
944  ntc.mtype[i] = vecPlay.mtype;
945  ntc.vecplay_ix[i] = vecPlay.ix;
946 #endif
947  Memb_list* ml = nt._ml_list[vecPlay.mtype];
948 
949  const std::vector<int>& array_dims = corenrn.get_array_dims()[vecPlay.mtype];
950 
951  auto padded_nodecount = nrn_soa_padded_size(ml->nodecount, Layout::SoA);
952  auto soaos_index = legacy2soaos_index(vecPlay.ix, array_dims);
953  vecPlay.ix = soaos2cnrn_index(soaos_index, array_dims, padded_nodecount, ml->_permute);
954 
955  nt._vecplay[i] = new VecPlayContinuous(ml->data + vecPlay.ix,
956  std::move(vecPlay.yvec),
957  std::move(vecPlay.tvec),
958  nullptr,
959  nt.id);
960  }
961 }
962 
963 void Phase2::populate(NrnThread& nt, const UserParams& userParams) {
965  ntc.file_id = userParams.gidgroups[nt.id];
966 
967  nt.ncell = n_real_cell;
968  nt.end = n_node;
970 
971 #if CHKPNTDEBUG
972  ntc.n_outputgids = n_output;
973  ntc.nmech = n_mech;
974 #endif
975 
976  /// Checkpoint in coreneuron is defined for both phase 1 and phase 2 since they are written
977  /// together
978  nt._ml_list = (Memb_list**) ecalloc_align(corenrn.get_memb_funcs().size(), sizeof(Memb_list*));
979 
980  auto& memb_func = corenrn.get_memb_funcs();
981 #if CHKPNTDEBUG
982  ntc.mlmap = new Memb_list_chkpnt*[memb_func.size()];
983  for (int i = 0; i < memb_func.size(); ++i) {
984  ntc.mlmap[i] = nullptr;
985  }
986 #endif
987 
988  nt.stream_id = 0;
989  nt.compute_gpu = 0;
992 
993 /* read_phase2 is being called from openmp region
994  * and hence we can set the stream equal to current thread id.
995  * In fact we could set gid as stream_id when we will have nrn threads
996  * greater than number of omp threads.
997  */
998 #if defined(_OPENMP)
999  nt.stream_id = omp_get_thread_num();
1000 #endif
1001 
1002  int shadow_rhs_cnt = 0;
1003  nt.shadow_rhs_cnt = 0;
1004 
1005  NrnThreadMembList* tml_last = nullptr;
1006  for (int i = 0; i < n_mech; ++i) {
1007  auto tml =
1008  create_tml(nt, i, memb_func[mech_types[i]], shadow_rhs_cnt, mech_types, nodecounts);
1009 
1010  nt._ml_list[tml->index] = tml->ml;
1011 
1012 #if CHKPNTDEBUG
1013  Memb_list_chkpnt* mlc = new Memb_list_chkpnt;
1014  ntc.mlmap[tml->index] = mlc;
1015 #endif
1016 
1017  if (nt.tml) {
1018  tml_last->next = tml;
1019  } else {
1020  nt.tml = tml;
1021  }
1022  tml_last = tml;
1023  }
1024 
1025  if (shadow_rhs_cnt) {
1026  nt._shadow_rhs = (double*) ecalloc_align(nrn_soa_padded_size(shadow_rhs_cnt, 0),
1027  sizeof(double));
1028  nt._shadow_d = (double*) ecalloc_align(nrn_soa_padded_size(shadow_rhs_cnt, 0),
1029  sizeof(double));
1030  nt.shadow_rhs_cnt = shadow_rhs_cnt;
1031  }
1032 
1033  nt.mapping = nullptr; // section segment mapping
1034 
1035  nt._nidata = n_idata;
1036  if (nt._nidata)
1037  nt._idata = (int*) ecalloc(nt._nidata, sizeof(int));
1038  else
1039  nt._idata = nullptr;
1040  // see patternstim.cpp
1041  int extra_nv = (&nt == nrn_threads) ? nrn_extra_thread0_vdata : 0;
1042  nt._nvdata = n_vdata;
1043  if (nt._nvdata + extra_nv)
1044  nt._vdata = (void**) ecalloc_align(nt._nvdata + extra_nv, sizeof(void*));
1045  else
1046  nt._vdata = nullptr;
1047 
1048  // The data format begins with the matrix data
1049  int n_data_padded = nrn_soa_padded_size(nt.end, SOA_LAYOUT);
1050  nt._data = _data;
1051  nt._actual_rhs = nt._data + 0 * n_data_padded;
1052  nt._actual_d = nt._data + 1 * n_data_padded;
1053  nt._actual_a = nt._data + 2 * n_data_padded;
1054  nt._actual_b = nt._data + 3 * n_data_padded;
1055  nt._actual_v = nt._data + 4 * n_data_padded;
1056  nt._actual_area = nt._data + 5 * n_data_padded;
1057  nt._actual_diam = n_diam ? nt._data + 6 * n_data_padded : nullptr;
1058 
1059  size_t offset = 6 * n_data_padded;
1060  if (n_diam) {
1061  // in the rare case that a mechanism has dparam with diam semantics
1062  // then actual_diam array added after matrix in nt._data
1063  // Generally wasteful since only a few diam are pointed to.
1064  // Probably better to move the diam semantics to the p array of the mechanism
1065  offset += n_data_padded;
1066  }
1067 
1068  // Memb_list.data points into the nt._data array.
1069  // Also count the number of Point_process
1070  int num_point_process = 0;
1071  for (auto tml = nt.tml; tml; tml = tml->next) {
1072  Memb_list* ml = tml->ml;
1073  int type = tml->index;
1074  int layout = corenrn.get_mech_data_layout()[type];
1075  int n = ml->nodecount;
1076  int sz = nrn_prop_param_size_[type];
1077  offset = nrn_soa_byte_align(offset);
1078  ml->data = nt._data + offset;
1079  offset += nrn_soa_padded_size(n, layout) * sz;
1080  if (corenrn.get_pnt_map()[type] > 0) {
1081  num_point_process += n;
1082  }
1083  }
1084  nt.pntprocs = new Point_process[num_point_process]{}; // includes acell with and without gid
1086  nt._ndata = offset;
1087 
1088 
1089  // matrix info
1091 
1092 #if CHKPNTDEBUG
1093  ntc.parent = new int[nt.end];
1094  memcpy(ntc.parent, nt._v_parent_index, nt.end * sizeof(int));
1095  ntc.area = new double[nt.end];
1096  memcpy(ntc.area, nt._actual_area, nt.end * sizeof(double));
1097 #endif
1098 
1099  int synoffset = 0;
1100  std::vector<int> pnt_offset(memb_func.size());
1101 
1102  // All the mechanism data and pdata.
1103  // Also fill in the pnt_offset
1104  // Complete spec of Point_process except for the acell presyn_ field.
1105  int itml = 0;
1106  for (auto tml = nt.tml; tml; tml = tml->next, ++itml) {
1107  int type = tml->index;
1108  Memb_list* ml = tml->ml;
1109  int n = ml->nodecount;
1110  int szp = nrn_prop_param_size_[type];
1111  int szdp = nrn_prop_dparam_size_[type];
1112 
1113  const std::vector<int>& array_dims = corenrn.get_array_dims()[type];
1114  int layout = corenrn.get_mech_data_layout()[type];
1115 
1116  ml->nodeindices = (int*) ecalloc_align(ml->nodecount, sizeof(int));
1117  std::copy(tmls[itml].nodeindices.begin(), tmls[itml].nodeindices.end(), ml->nodeindices);
1118 
1119  mech_data_layout_transform<double>(ml->data, n, array_dims, layout);
1120 
1121  if (szdp) {
1122  ml->pdata = (int*) ecalloc_align(nrn_soa_padded_size(n, layout) * szdp, sizeof(int));
1123  std::copy(tmls[itml].pdata.begin(), tmls[itml].pdata.end(), ml->pdata);
1124  mech_data_layout_transform<int>(ml->pdata, n, szdp, layout);
1125 
1126 #if CHKPNTDEBUG // Not substantive. Only for debugging.
1127  Memb_list_chkpnt* mlc = ntc.mlmap[type];
1128  mlc->pdata_not_permuted = (int*) coreneuron::ecalloc_align(n * szdp, sizeof(int));
1129  if (layout == Layout::AoS) { // only copy
1130  for (int i = 0; i < n; ++i) {
1131  for (int j = 0; j < szdp; ++j) {
1132  mlc->pdata_not_permuted[i * szdp + j] = ml->pdata[i * szdp + j];
1133  }
1134  }
1135  } else if (layout == Layout::SoA) { // transpose and unpad
1136  int align_cnt = nrn_soa_padded_size(n, layout);
1137  for (int i = 0; i < n; ++i) {
1138  for (int j = 0; j < szdp; ++j) {
1139  mlc->pdata_not_permuted[i * szdp + j] = ml->pdata[i + j * align_cnt];
1140  }
1141  }
1142  }
1143 #endif
1144  } else {
1145  ml->pdata = nullptr;
1146  }
1147  if (corenrn.get_pnt_map()[type] > 0) { // POINT_PROCESS mechanism including acell
1148  int cnt = ml->nodecount;
1149  Point_process* pnt = nullptr;
1150  pnt = nt.pntprocs + synoffset;
1151  pnt_offset[type] = synoffset;
1152  synoffset += cnt;
1153  for (int i = 0; i < cnt; ++i) {
1154  Point_process* pp = pnt + i;
1155  pp->_type = type;
1156  pp->_i_instance = i;
1157  nt._vdata[ml->pdata[nrn_i_layout(i, cnt, 1, szdp, layout)]] = pp;
1158  pp->_tid = nt.id;
1159  }
1160  }
1161 
1162  auto& r = tmls[itml].nmodlrandom;
1163  if (r.size()) {
1164  size_t ix{};
1165  uint32_t n_randomvar = r[ix++];
1166  assert(r.size() == 1 + n_randomvar + 5 * n_randomvar * n);
1167  std::vector<uint32_t> indices(n_randomvar);
1168  for (uint32_t i = 0; i < n_randomvar; ++i) {
1169  indices[i] = r[ix++];
1170  }
1171  int cnt = ml->nodecount;
1172  for (auto index: indices) {
1173  // should we also verify that index on corenrn side same as on nrn side?
1174  // sonarcloud thinks ml_pdata can be nullptr, so ...
1175  assert(index >= 0 && index < szdp);
1176  for (int i = 0; i < n; ++i) {
1177  nrnran123_State* state = nrnran123_newstream3(r[ix], r[ix + 1], r[ix + 2]);
1178  nrnran123_setseq(state, r[ix + 3], char(r[ix + 4]));
1179  ix += 5;
1180  int ipd = ml->pdata[nrn_i_layout(i, cnt, index, szdp, layout)];
1181  assert(ipd >= 0 && ipd < n_vdata + extra_nv);
1182  nt._vdata[ipd] = state;
1183  }
1184  }
1185  }
1186  }
1187 
1188  // pnt_offset needed for SelfEvent transfer from NEURON. Not needed on GPU.
1189  // Ugh. Related but not same as NetReceiveBuffer._pnt_offset
1190  nt._pnt_offset = pnt_offset;
1191 
1193 
1194  /* if desired, apply the node permutation. This involves permuting
1195  at least the node parameter arrays for a, b, and area (and diam) and all
1196  integer vector values that index into nodes. This could have been done
1197  when originally filling the arrays with AoS ordered data, but can also
1198  be done now, after the SoA transformation. The latter has the advantage
1199  that the present order is consistent with all the layout values. Note
1200  that after this portion of the permutation, a number of other node index
1201  vectors will be read and will need to be permuted as well in subsequent
1202  sections of this function.
1203  */
1205  nt._permute = interleave_order(nt.id, nt.ncell, nt.end, nt._v_parent_index);
1206  }
1207  if (nt._permute) {
1208  int* p = nt._permute;
1209  permute_data(nt._actual_a, nt.end, p);
1210  permute_data(nt._actual_b, nt.end, p);
1211  permute_data(nt._actual_area, nt.end, p);
1213  nt.end,
1214  p); // need if restore or finitialize does not initialize voltage
1215  if (nt._actual_diam) {
1216  permute_data(nt._actual_diam, nt.end, p);
1217  }
1218  // index values change as well as ordering
1219  permute_ptr(nt._v_parent_index, nt.end, p);
1220  node_permute(nt._v_parent_index, nt.end, p);
1221 
1222 #if CORENRN_DEBUG
1223  for (int i = 0; i < nt.end; ++i) {
1224  printf("parent[%d] = %d\n", i, nt._v_parent_index[i]);
1225  }
1226 #endif
1227 
1228  // specify the ml->_permute and sort the nodeindices
1229  // Have to calculate all the permute before updating pdata in case
1230  // POINTER to data of other mechanisms exist.
1231  for (auto tml = nt.tml; tml; tml = tml->next) {
1232  if (tml->ml->nodeindices) { // not artificial
1233  permute_nodeindices(tml->ml, p);
1234  }
1235  }
1236  for (auto tml = nt.tml; tml; tml = tml->next) {
1237  if (tml->ml->nodeindices) { // not artificial
1238  permute_ml(tml->ml, tml->index, nt);
1239  }
1240  }
1241 
1242  // permute the Point_process._i_instance
1243  for (int i = 0; i < nt.n_pntproc; ++i) {
1244  Point_process& pp = nt.pntprocs[i];
1245  Memb_list* ml = nt._ml_list[pp._type];
1246  if (ml->_permute) {
1247  pp._i_instance = ml->_permute[pp._i_instance];
1248  }
1249  }
1250  }
1251 
1253 
1255 
1256  // for fast watch statement checking
1257  // setup a list of types that have WATCH statement
1258  {
1259  int sz = 0; // count the types with WATCH
1260  for (auto tml = nt.tml; tml; tml = tml->next) {
1261  if (corenrn.get_watch_check()[tml->index]) {
1262  ++sz;
1263  }
1264  }
1265  if (sz) {
1266  nt._watch_types = (int*) ecalloc(sz + 1, sizeof(int)); // nullptr terminated
1267  sz = 0;
1268  for (auto tml = nt.tml; tml; tml = tml->next) {
1269  if (corenrn.get_watch_check()[tml->index]) {
1270  nt._watch_types[sz++] = tml->index;
1271  }
1272  }
1273  }
1274  }
1275  auto& pnttype2presyn = corenrn.get_pnttype2presyn();
1277  // create the nt.pnt2presyn_ix array of arrays.
1278  nt.pnt2presyn_ix = (int**) ecalloc(nrn_has_net_event_.size(), sizeof(int*));
1279  for (size_t i = 0; i < nrn_has_net_event_.size(); ++i) {
1281  if (ml && ml->nodecount > 0) {
1282  nt.pnt2presyn_ix[i] = (int*) ecalloc(ml->nodecount, sizeof(int));
1283  }
1284  }
1285 
1286  // Real cells are at the beginning of the nt.presyns followed by
1287  // acells (with and without gids mixed together)
1288  // Here we associate the real cells with voltage pointers and
1289  // acell PreSyn with the Point_process.
1290  // nt.presyns order same as output_vindex order
1291 #if CHKPNTDEBUG
1292  ntc.output_vindex = new int[nt.n_presyn];
1293  memcpy(ntc.output_vindex, output_vindex.data(), nt.n_presyn * sizeof(int));
1294 #endif
1295  if (nt._permute) {
1296  // only indices >= 0 (i.e. _actual_v indices) will be changed.
1297  node_permute(output_vindex.data(), nt.n_presyn, nt._permute);
1298  }
1299 #if CHKPNTDEBUG
1300  ntc.output_threshold = new double[n_real_output];
1301  memcpy(ntc.output_threshold, output_threshold.data(), n_real_output * sizeof(double));
1302 #endif
1303 
1304  for (int i = 0; i < nt.n_presyn; ++i) { // real cells
1305  PreSyn* ps = nt.presyns + i;
1306 
1307  int ix = output_vindex[i];
1308  if (ix == -1 && i < n_real_output) { // real cell without a presyn
1309  continue;
1310  }
1311  if (ix < 0) {
1312  ix = -ix;
1313  int index = ix / 1000;
1314  int type = ix % 1000;
1315  Point_process* pnt = nt.pntprocs + (pnt_offset[type] + index);
1316  ps->pntsrc_ = pnt;
1317  // pnt->_presyn = ps;
1318  int ip2ps = pnttype2presyn[pnt->_type];
1319  if (ip2ps >= 0) {
1320  nt.pnt2presyn_ix[ip2ps][pnt->_i_instance] = i;
1321  }
1322  if (ps->gid_ < 0) {
1323  ps->gid_ = -1;
1324  }
1325  } else {
1326  assert(ps->gid_ > -1);
1327  ps->thvar_index_ = ix; // index into _actual_v
1328  assert(ix < nt.end);
1329  ps->threshold_ = output_threshold[i];
1330  }
1331  }
1332 
1333  // initial net_send_buffer size about 1% of number of presyns
1334  // nt._net_send_buffer_size = nt.ncell/100 + 1;
1335  // but, to avoid reallocation complexity on GPU ...
1337  nt._net_send_buffer = (int*) ecalloc_align(nt._net_send_buffer_size, sizeof(int));
1338 
1339  int nnetcon = nt.n_netcon;
1340 
1341  // it may happen that Point_process structures will be made unnecessary
1342  // by factoring into NetCon.
1343 
1344 #if CHKPNTDEBUG
1345  ntc.pnttype = new int[nnetcon];
1346  ntc.pntindex = new int[nnetcon];
1347  memcpy(ntc.pnttype, pnttype.data(), nnetcon * sizeof(int));
1348  memcpy(ntc.pntindex, pntindex.data(), nnetcon * sizeof(int));
1349 #endif
1350  for (int i = 0; i < nnetcon; ++i) {
1351  int type = pnttype[i];
1352  if (type > 0) {
1353  int index = pnt_offset[type] + pntindex[i]; /// Potentially uninitialized pnt_offset[],
1354  /// check for previous assignments
1355  NetCon& nc = nt.netcons[i];
1356  nc.target_ = nt.pntprocs + index;
1357  nc.active_ = true;
1358  }
1359  }
1360 
1361  handle_weights(nt, nnetcon, ntc);
1362 
1363  get_info_from_bbcore(nt, memb_func, ntc);
1364 
1365  set_vec_play(nt, ntc);
1366 
1367  if (!events.empty()) {
1368  userParams.checkPoints.restore_tqueue(nt, *this);
1369  }
1370 
1371  set_net_send_buffer(nt._ml_list, pnt_offset);
1372 }
1373 } // namespace coreneuron
int * nrn_has_net_event_
Definition: init.cpp:161
Definition: utility.h:81
void restore_tqueue(NrnThread &, const Phase2 &p2)
auto & get_net_buf_send_type()
Definition: coreneuron.hpp:155
auto & get_different_mechanism_type()
Definition: coreneuron.hpp:139
T * read_array(T *p, size_t count)
Read an integer array of fixed length.
std::vector< T > read_vector(size_t count)
void record_checkpoint()
Record current chkpnt state.
bool eof()
nothing more to read
int read_int()
Parse a single integer entry.
union coreneuron::NetCon::@0 u
Point_process * target_
Definition: netcon.hpp:49
void set_net_send_buffer(Memb_list **ml_list, const std::vector< int > &pnt_offset)
Definition: phase2.cpp:511
std::vector< int > pntindex
Definition: phase2.hpp:125
std::vector< int > mech_types
Definition: phase2.hpp:99
std::vector< int > nodecounts
Definition: phase2.hpp:100
int * v_parent_index
Definition: phase2.hpp:103
std::vector< double > weights
Definition: phase2.hpp:126
void pdata_relocation(const NrnThread &nt, const std::vector< Memb_func > &memb_func)
Definition: phase2.cpp:639
void set_vec_play(NrnThread &nt, NrnThreadChkpnt &ntc)
Definition: phase2.cpp:925
void handle_weights(NrnThread &nt, int n_netcon, NrnThreadChkpnt &ntc)
Definition: phase2.cpp:835
std::vector< TML > tmls
Definition: phase2.hpp:121
double * _data
Definition: phase2.hpp:111
void get_info_from_bbcore(NrnThread &nt, const std::vector< Memb_func > &memb_func, NrnThreadChkpnt &ntc)
Definition: phase2.cpp:867
void check_mechanism()
Check if MOD file used between NEURON and CoreNEURON is same.
Definition: phase2.cpp:470
void transform_int_data(int elem0, int nodecount, int *pdata, int i, int dparam_size, int layout, int n_node_)
Perform in memory transformation between AoS<>SoA for integer data.
Definition: phase2.cpp:496
void restore_events(FileHandler &F)
Definition: phase2.cpp:548
std::vector< double > output_threshold
Definition: phase2.hpp:123
std::vector< int > preSynConditionEventFlags
Definition: phase2.hpp:31
std::vector< int > pnttype
Definition: phase2.hpp:124
std::vector< VecPlayContinuous_ > vec_play_continuous
Definition: phase2.hpp:68
void read_direct(int thread_id, const NrnThread &nt)
Definition: phase2.cpp:279
std::vector< int > output_vindex
Definition: phase2.hpp:122
void read_file(FileHandler &F, const NrnThread &nt)
Definition: phase2.cpp:126
void set_dependencies(const NrnThread &nt, const std::vector< Memb_func > &memb_func)
Definition: phase2.cpp:765
std::vector< double > delay
Definition: phase2.hpp:127
std::vector< std::pair< int, std::shared_ptr< EventTypeBase > > > events
Definition: phase2.hpp:71
void populate(NrnThread &nt, const UserParams &userParams)
Definition: phase2.cpp:963
void fill_before_after_lists(NrnThread &nt, const std::vector< Memb_func > &memb_func)
Definition: phase2.cpp:606
Point_process * pntsrc_
Definition: netcon.hpp:113
#define cnt
Definition: tqueue.hpp:44
#define nodeindices
Definition: md1redef.h:35
#define area
Definition: md1redef.h:12
#define v
Definition: md1redef.h:11
#define data
Definition: md1redef.h:36
#define nodecount
Definition: md1redef.h:39
#define weights
Definition: md1redef.h:42
#define i
Definition: md1redef.h:19
#define pdata
Definition: md1redef.h:37
#define VecPlayContinuousType
Definition: vrecitem.h:17
#define PlayRecordEventType
Definition: vrecitem.h:18
#define SOA_LAYOUT
Definition: data_layout.hpp:11
static RNG::key_type k
Definition: nrnran123.cpp:9
void nrnran123_setseq(nrnran123_State *s, std::uint32_t seq, char which)
Set a Random123 sequence for a sequnece ID and which selector.
Definition: nrnran123.cpp:55
#define assert(ex)
Definition: hocassrt.h:24
bool nrn_semantics_is_ion(int i)
Definition: ion_semantics.h:6
int nrn_semantics_ion_type(int i)
Definition: ion_semantics.h:12
#define BEFORE_AFTER_SIZE
Definition: membfunc.hpp:72
printf
Definition: extdef.h:5
void move(Item *q1, Item *q2, Item *q3)
Definition: list.cpp:200
THIS FILE IS AUTO GENERATED DONT MODIFY IT.
NrnThread * nrn_threads
Definition: multicore.cpp:56
void setup_fornetcon_info(NrnThread &nt)
If FOR_NETCON in use, setup NrnThread fornetcon related info.
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.
int nrn_i_layout(int icnt, int cnt, int isz, int sz, int layout)
This function return the index in a flat array of a matrix coordinate (icnt, isz).
fixed_vector< double > IvocVect
Definition: ivocvect.hpp:72
nrnran123_State * nrnran123_newstream3(uint32_t id1, uint32_t id2, uint32_t id3, bool use_unified_memory)
Allocate a new Random123 stream.
Definition: nrnran123.cpp:172
void nrn_abort(int errcode)
Definition: utils.cpp:13
void * ecalloc_align(size_t n, size_t size, size_t alignment)
int nrn_extra_thread0_vdata
Definition: patternstim.cpp:46
int Datum
Definition: nrnconf.h:23
int interleave_permute_type
const char * nrn_get_mechname(int type)
Definition: mk_mech.cpp:152
void mech_data_layout_transform(T *data, int cnt, const std::vector< int > &array_dims, int layout)
Definition: phase2.cpp:91
int nrn_mech_depend(int type, int *dependencies)
std::array< int, 3 > legacy2soaos_index(int legacy_index, const std::vector< int > &array_dims)
Split a legacy index into the three SoAoS indices.
CoreNeuron corenrn
Definition: multicore.cpp:53
void * ecalloc(size_t n, size_t size)
Definition: nrnoc_aux.cpp:85
static void * emalloc(size_t size)
Definition: mpispike.cpp:30
NrnThreadChkpnt * nrnthread_chkpnt
size_t nrn_soa_byte_align(size_t size)
return the new offset considering the byte aligment settings
int nrn_soa_padded_size(int cnt, int layout)
calculate size after padding for specific memory layout
NrnThreadMembList * create_tml(NrnThread &nt, int mech_id, Memb_func &memb_func, int &shadow_rhs_cnt, const std::vector< int > &mech_types, const std::vector< int > &nodecounts)
Definition: multicore.cpp:64
void permute_nodeindices(Memb_list *ml, int *permute)
void permute_ml(Memb_list *ml, int type, NrnThread &nt)
void permute_data(double *vec, int n, int *p)
std::vector< int > interleave_order(int ith, int ncell, int nnode, int *parent)
Function that performs the permutation of the cells such that the execution threads access coalesced ...
Definition: cellorder.cpp:348
int ii
Definition: cellorder.cpp:631
void permute_ptr(int *vec, int n, int *p)
#define NetParEventType
Definition: netcon.hpp:27
#define PreSynType
Definition: netcon.hpp:26
#define SelfEventType
Definition: netcon.hpp:25
#define NetConType
Definition: netcon.hpp:24
#define nrn_assert(x)
assert()-like macro, independent of NDEBUG status
Definition: nrn_assert.h:33
int * nrn_prop_dparam_size_
Definition: init.cpp:163
int const size_t const size_t n
Definition: nrngsl.h:10
size_t p
size_t j
s
Definition: multisend.cpp:521
int * nrn_prop_param_size_
Definition: init.cpp:162
short index
Definition: cabvars.h:11
std::vector< Memb_func > memb_func
Definition: init.cpp:145
short type
Definition: cabvars.h:10
static int pcnt
Definition: axis.cpp:425
int(* nrn2core_get_dat2_mech_)(int tid, size_t i, int dsz_inst, int *&nodeindices, double *&data, int *&pdata, std::vector< uint32_t > &nmodlrandom, std::vector< int > &pointer2type)
Definition: phase2.cpp:48
int(* nrn2core_get_dat2_3_)(int tid, int nweight, int *&output_vindex, double *&output_threshold, int *&netcon_pnttype, int *&netcon_pntindex, double *&weights, double *&delays)
Definition: phase2.cpp:57
int(* nrn2core_get_dat2_corepointer_mech_)(int tid, int type, int &icnt, int &dcnt, int *&iarray, double *&darray)
Definition: phase2.cpp:68
int(* nrn2core_get_dat2_vecplay_inst_)(int tid, int i, int &vptype, int &mtype, int &ix, int &sz, double *&yvec, double *&tvec, int &last_index, int &discon_index, int &ubound_index)
Definition: phase2.cpp:77
int(* nrn2core_get_dat2_corepointer_)(int tid, int &n)
Definition: phase2.cpp:66
int(* nrn2core_get_dat2_2_)(int tid, int *&v_parent_index, double *&a, double *&b, double *&area, double *&v, double *&diamvec)
Definition: phase2.cpp:40
int(* nrn2core_get_dat2_1_)(int tid, int &n_real_cell, int &ngid, int &n_real_gid, int &nnode, int &ndiam, int &nmech, int *&tml_index, int *&ml_nodecount, int &nidata, int &nvdata, int &nweight)
Definition: phase2.cpp:27
int(* nrn2core_get_dat2_vecplay_)(int tid, std::vector< int > &indices)
Definition: phase2.cpp:75
ThreadDatum * _thread
Definition: mechanism.hpp:141
NetSendBuffer_t * _net_send_buffer
Definition: mechanism.hpp:143
NetReceiveBuffer_t * _net_receive_buffer
Definition: mechanism.hpp:142
NrnThreadBAList * next
Definition: multicore.hpp:49
NrnThreadBAList * tbl[BEFORE_AFTER_SIZE]
Definition: multicore.hpp:133
std::vector< int > _pnt_offset
Definition: multicore.hpp:154
Memb_list ** _ml_list
Definition: multicore.hpp:81
NrnThreadMembList * tml
Definition: multicore.hpp:80
Point_process * pntprocs
Definition: multicore.hpp:82
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 *const gidgroups
Array of cell group numbers (indices)
Definition: user_params.hpp:35
CheckPoints & checkPoints
Definition: user_params.hpp:41