NEURON
mpispike.cpp
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 #include <stdio.h>
3 #include <stdlib.h>
4 #include <stddef.h>
5 #include <assert.h>
6 
7 /* do not want the redef in the dynamic load case */
8 #include <nrnmpiuse.h>
9 
10 #if NRNMPI_DYNAMICLOAD
11 #include <nrnmpi_dynam.h>
12 #endif
13 
14 #include <nrnmpi.h>
15 #include <hocdec.h>
16 
17 #if NRNMPI
18 #include "nrnmpidec.h"
19 #include "nrnmpi_impl.h"
20 #include "mpispike.h"
21 #include <mpi.h>
22 
23 extern void nrnbbs_context_wait();
24 
25 static int np;
26 static int* displs;
27 static int* byteovfl; /* for the compressed transfer method */
28 static MPI_Datatype spike_type;
29 
30 static void pgvts_op(double* in, double* inout, int* len, MPI_Datatype* dptr);
31 static MPI_Op mpi_pgvts_op;
32 
33 static void make_spike_type() {
34  NRNMPI_Spike s;
35  int block_lengths[2];
36  MPI_Aint displacements[2];
37  MPI_Aint addresses[3];
38  MPI_Datatype typelist[2];
39 
40  typelist[0] = MPI_INT;
41  typelist[1] = MPI_DOUBLE;
42 
43  block_lengths[0] = block_lengths[1] = 1;
44 
45  MPI_Get_address(&s, &addresses[0]);
46  MPI_Get_address(&(s.gid), &addresses[1]);
47  MPI_Get_address(&(s.spiketime), &addresses[2]);
48 
49  displacements[0] = addresses[1] - addresses[0];
50  displacements[1] = addresses[2] - addresses[0];
51 
52  MPI_Type_create_struct(2, block_lengths, displacements, typelist, &spike_type);
53  MPI_Type_commit(&spike_type);
54 
55  MPI_Op_create((MPI_User_function*) pgvts_op, 1, &mpi_pgvts_op);
56 }
57 
58 void nrnmpi_spike_initialize() {
59  make_spike_type();
60 }
61 
62 #if nrn_spikebuf_size > 0
63 
64 static MPI_Datatype spikebuf_type;
65 
66 static void make_spikebuf_type() {
67  NRNMPI_Spikebuf s;
68  int block_lengths[3];
69  MPI_Aint displacements[3];
70  MPI_Aint addresses[4];
71  MPI_Datatype typelist[3];
72 
73  typelist[0] = MPI_INT;
74  typelist[1] = MPI_INT;
75  typelist[2] = MPI_DOUBLE;
76 
77  block_lengths[0] = 1;
78  block_lengths[1] = nrn_spikebuf_size;
79  block_lengths[2] = nrn_spikebuf_size;
80 
81  MPI_Get_address(&s, &addresses[0]);
82  MPI_Get_address(&(s.nspike), &addresses[1]);
83  MPI_Get_address(&(s.gid[0]), &addresses[2]);
84  MPI_Get_address(&(s.spiketime[0]), &addresses[3]);
85 
86  displacements[0] = addresses[1] - addresses[0];
87  displacements[1] = addresses[2] - addresses[0];
88  displacements[2] = addresses[3] - addresses[0];
89 
90  MPI_Type_create_struct(3, block_lengths, displacements, typelist, &spikebuf_type);
91  MPI_Type_commit(&spikebuf_type);
92 }
93 #endif
94 
95 int nrnmpi_spike_exchange() {
96  int i, n, novfl, n1;
97  if (!displs) {
98  np = nrnmpi_numprocs;
99  displs = (int*) hoc_Emalloc(np * sizeof(int));
100  hoc_malchk();
101  displs[0] = 0;
102 #if nrn_spikebuf_size > 0
103  make_spikebuf_type();
104 #endif
105  }
107 #if nrn_spikebuf_size == 0
108  MPI_Allgather(&nout_, 1, MPI_INT, nin_, 1, MPI_INT, nrnmpi_comm);
109  n = nin_[0];
110  for (i = 1; i < np; ++i) {
111  displs[i] = n;
112  n += nin_[i];
113  }
114  if (n) {
115  if (icapacity_ < n) {
116  icapacity_ = n + 10;
117  free(spikein_);
119  hoc_malchk();
120  }
121  MPI_Allgatherv(
122  spikeout_, nout_, spike_type, spikein_, nin_, displs, spike_type, nrnmpi_comm);
123  }
124 #else
125  MPI_Allgather(spbufout_, 1, spikebuf_type, spbufin_, 1, spikebuf_type, nrnmpi_comm);
126  novfl = 0;
127  n = spbufin_[0].nspike;
128  if (n > nrn_spikebuf_size) {
129  nin_[0] = n - nrn_spikebuf_size;
130  novfl += nin_[0];
131  } else {
132  nin_[0] = 0;
133  }
134  for (i = 1; i < np; ++i) {
135  displs[i] = novfl;
136  n1 = spbufin_[i].nspike;
137  n += n1;
138  if (n1 > nrn_spikebuf_size) {
139  nin_[i] = n1 - nrn_spikebuf_size;
140  novfl += nin_[i];
141  } else {
142  nin_[i] = 0;
143  }
144  }
145  if (novfl) {
146  if (icapacity_ < novfl) {
147  icapacity_ = novfl + 10;
148  free(spikein_);
150  hoc_malchk();
151  }
153  MPI_Allgatherv(spikeout_, n1, spike_type, spikein_, nin_, displs, spike_type, nrnmpi_comm);
154  }
155  ovfl_ = novfl;
156 #endif
157  return n;
158 }
159 
160 /*
161 The compressed spike format is restricted to the fixed step method and is
162 a sequence of unsigned char.
163 nspike = buf[0]*256 + buf[1]
164 a sequence of spiketime, localgid pairs. There are nspike of them.
165  spiketime is relative to the last transfer time in units of dt.
166  note that this requires a mindelay < 256*dt.
167  localgid is an unsigned int, unsigned short,
168  or unsigned char in size depending on the range and thus takes
169  4, 2, or 1 byte respectively. To be machine independent we do our
170  own byte coding. When the localgid range is smaller than the true
171  gid range, the gid->PreSyn are remapped into
172  hostid specific maps. If there are not many holes, i.e just about every
173  spike from a source machine is delivered to some cell on a
174  target machine, then instead of a hash map, a vector is used.
175 The allgather sends the first part of the buf and the allgatherv buffer
176 sends any overflow.
177 */
178 int nrnmpi_spike_exchange_compressed() {
179  int i, novfl, n, ntot, idx, bs, bstot; /* n is #spikes, bs is #byte overflow */
180  if (!displs) {
181  np = nrnmpi_numprocs;
182  displs = (int*) hoc_Emalloc(np * sizeof(int));
183  hoc_malchk();
184  displs[0] = 0;
185  }
186  if (!byteovfl) {
187  byteovfl = (int*) hoc_Emalloc(np * sizeof(int));
188  hoc_malchk();
189  }
191 
192  MPI_Allgather(
193  spfixout_, ag_send_size_, MPI_BYTE, spfixin_, ag_send_size_, MPI_BYTE, nrnmpi_comm);
194  novfl = 0;
195  ntot = 0;
196  bstot = 0;
197  for (i = 0; i < np; ++i) {
198  displs[i] = bstot;
199  idx = i * ag_send_size_;
200  n = spfixin_[idx++] * 256;
201  n += spfixin_[idx++];
202  ntot += n;
203  nin_[i] = n;
204  if (n > ag_send_nspike_) {
205  bs = 2 + n * (1 + localgid_size_) - ag_send_size_;
206  byteovfl[i] = bs;
207  bstot += bs;
208  novfl += n - ag_send_nspike_;
209  } else {
210  byteovfl[i] = 0;
211  }
212  }
213  if (novfl) {
214  if (ovfl_capacity_ < novfl) {
215  ovfl_capacity_ = novfl + 10;
216  free(spfixin_ovfl_);
217  spfixin_ovfl_ = (unsigned char*) hoc_Emalloc(ovfl_capacity_ * (1 + localgid_size_) *
218  sizeof(unsigned char));
219  hoc_malchk();
220  }
221  bs = byteovfl[nrnmpi_myid];
222  /*
223  note that the spfixout_ buffer is one since the overflow
224  is contiguous to the first part. But the spfixin_ovfl_ is
225  completely separate from the spfixin_ since the latter
226  dynamically changes its size during a run.
227  */
228  MPI_Allgatherv(spfixout_ + ag_send_size_,
229  bs,
230  MPI_BYTE,
232  byteovfl,
233  displs,
234  MPI_BYTE,
235  nrnmpi_comm);
236  }
237  ovfl_ = novfl;
238  return ntot;
239 }
240 
241 double nrnmpi_mindelay(double m) {
242  double result;
243  if (!nrnmpi_use) {
244  return m;
245  }
247  MPI_Allreduce(&m, &result, 1, MPI_DOUBLE, MPI_MIN, nrnmpi_comm);
248  return result;
249 }
250 
251 int nrnmpi_int_allmax(int x) {
252  int result;
253  if (nrnmpi_numprocs < 2) {
254  return x;
255  }
257  MPI_Allreduce(&x, &result, 1, MPI_INT, MPI_MAX, nrnmpi_comm);
258  return result;
259 }
260 
261 #define ALLTOALLV_SPARSE_TAG 101980
262 
263 /* Code derived from MPI_Alltoallv_sparse in MP-Gadget: https://github.com/MP-Gadget */
264 
265 static int MPI_Alltoallv_sparse(void* sendbuf,
266  int* sendcnts,
267  int* sdispls,
268  MPI_Datatype sendtype,
269  void* recvbuf,
270  int* recvcnts,
271  int* rdispls,
272  MPI_Datatype recvtype,
273  MPI_Comm comm) {
274  int status;
275  int myrank;
276  int nranks;
277  status = MPI_Comm_rank(comm, &myrank);
278  assert(status == MPI_SUCCESS);
279  status = MPI_Comm_size(comm, &nranks);
280  assert(status == MPI_SUCCESS);
281 
282  int rankp;
283  for (rankp = 0; nranks > (1 << rankp); rankp++)
284  ;
285 
286  ptrdiff_t lb;
287  ptrdiff_t send_elsize;
288  ptrdiff_t recv_elsize;
289 
290  status = MPI_Type_get_extent(sendtype, &lb, &send_elsize);
291  assert(status == MPI_SUCCESS);
292  status = MPI_Type_get_extent(recvtype, &lb, &recv_elsize);
293  assert(status == MPI_SUCCESS);
294 
295  MPI_Request* requests = (MPI_Request*) hoc_Emalloc(nranks * 2 * sizeof(MPI_Request));
296  hoc_malchk();
297  assert(requests != NULL);
298 
299  int ngrp;
300  int n_requests;
301  n_requests = 0;
302  for (ngrp = 0; ngrp < (1 << rankp); ngrp++) {
303  int target = myrank ^ ngrp;
304 
305  if (target >= nranks)
306  continue;
307  if (recvcnts[target] == 0)
308  continue;
309  status = MPI_Irecv((static_cast<char*>(recvbuf)) + recv_elsize * rdispls[target],
310  recvcnts[target],
311  recvtype,
312  target,
313  ALLTOALLV_SPARSE_TAG,
314  comm,
315  &requests[n_requests++]);
316  assert(status == MPI_SUCCESS);
317  }
318 
319  status = MPI_Barrier(comm);
320  assert(status == MPI_SUCCESS);
321 
322  for (ngrp = 0; ngrp < (1 << rankp); ngrp++) {
323  int target = myrank ^ ngrp;
324  if (target >= nranks)
325  continue;
326  if (sendcnts[target] == 0)
327  continue;
328  status = MPI_Isend((static_cast<char*>(sendbuf)) + send_elsize * sdispls[target],
329  sendcnts[target],
330  sendtype,
331  target,
332  ALLTOALLV_SPARSE_TAG,
333  comm,
334  &requests[n_requests++]);
335  assert(status == MPI_SUCCESS);
336  }
337 
338  status = MPI_Waitall(n_requests, requests, MPI_STATUSES_IGNORE);
339  assert(status == MPI_SUCCESS);
340  free(requests);
341 
342  status = MPI_Barrier(comm);
343  assert(status == MPI_SUCCESS);
344 
345  return MPI_SUCCESS;
346 }
347 
348 
349 extern void nrnmpi_dbl_alltoallv_sparse(double* s,
350  int* scnt,
351  int* sdispl,
352  double* r,
353  int* rcnt,
354  int* rdispl) {
355  MPI_Alltoallv_sparse(s, scnt, sdispl, MPI_DOUBLE, r, rcnt, rdispl, MPI_DOUBLE, nrnmpi_comm);
356 }
357 extern void nrnmpi_int_alltoallv_sparse(int* s,
358  int* scnt,
359  int* sdispl,
360  int* r,
361  int* rcnt,
362  int* rdispl) {
363  MPI_Alltoallv_sparse(s, scnt, sdispl, MPI_INT, r, rcnt, rdispl, MPI_INT, nrnmpi_comm);
364 }
365 
366 extern void nrnmpi_long_alltoallv_sparse(int64_t* s,
367  int* scnt,
368  int* sdispl,
369  int64_t* r,
370  int* rcnt,
371  int* rdispl) {
372  MPI_Alltoallv_sparse(s, scnt, sdispl, MPI_INT64_T, r, rcnt, rdispl, MPI_INT64_T, nrnmpi_comm);
373 }
374 
375 
376 extern void nrnmpi_int_gather(int* s, int* r, int cnt, int root) {
377  MPI_Gather(s, cnt, MPI_INT, r, cnt, MPI_INT, root, nrnmpi_comm);
378 }
379 
380 extern void nrnmpi_int_gatherv(int* s, int scnt, int* r, int* rcnt, int* rdispl, int root) {
381  MPI_Gatherv(s, scnt, MPI_INT, r, rcnt, rdispl, MPI_INT, root, nrnmpi_comm);
382 }
383 
384 extern void nrnmpi_char_gatherv(char* s, int scnt, char* r, int* rcnt, int* rdispl, int root) {
385  MPI_Gatherv(s, scnt, MPI_CHAR, r, rcnt, rdispl, MPI_CHAR, root, nrnmpi_comm);
386 }
387 
388 extern void nrnmpi_int_scatter(int* s, int* r, int cnt, int root) {
389  MPI_Scatter(s, cnt, MPI_INT, r, cnt, MPI_INT, root, nrnmpi_comm);
390 }
391 
392 extern void nrnmpi_char_scatterv(char* s, int* scnt, int* sdispl, char* r, int rcnt, int root) {
393  MPI_Scatterv(s, scnt, sdispl, MPI_CHAR, r, rcnt, MPI_CHAR, root, nrnmpi_comm);
394 }
395 
396 extern void nrnmpi_int_alltoall(int* s, int* r, int n) {
397  MPI_Alltoall(s, n, MPI_INT, r, n, MPI_INT, nrnmpi_comm);
398 }
399 
400 extern void nrnmpi_int_alltoallv(int* s, int* scnt, int* sdispl, int* r, int* rcnt, int* rdispl) {
401  MPI_Alltoallv(s, scnt, sdispl, MPI_INT, r, rcnt, rdispl, MPI_INT, nrnmpi_comm);
402 }
403 
404 extern void nrnmpi_long_alltoallv(int64_t* s,
405  int* scnt,
406  int* sdispl,
407  int64_t* r,
408  int* rcnt,
409  int* rdispl) {
410  MPI_Alltoallv(s, scnt, sdispl, MPI_INT64_T, r, rcnt, rdispl, MPI_INT64_T, nrnmpi_comm);
411 }
412 
413 extern void nrnmpi_dbl_alltoallv(double* s,
414  int* scnt,
415  int* sdispl,
416  double* r,
417  int* rcnt,
418  int* rdispl) {
419  MPI_Alltoallv(s, scnt, sdispl, MPI_DOUBLE, r, rcnt, rdispl, MPI_DOUBLE, nrnmpi_comm);
420 }
421 
422 extern void nrnmpi_char_alltoallv(char* s,
423  int* scnt,
424  int* sdispl,
425  char* r,
426  int* rcnt,
427  int* rdispl) {
428  MPI_Alltoallv(s, scnt, sdispl, MPI_CHAR, r, rcnt, rdispl, MPI_CHAR, nrnmpi_comm);
429 }
430 
431 /* following are for the partrans */
432 
433 void nrnmpi_int_allgather(int* s, int* r, int n) {
434  MPI_Allgather(s, n, MPI_INT, r, n, MPI_INT, nrnmpi_comm);
435 }
436 
437 void nrnmpi_int_allgather_inplace(int* srcdest, int n) {
438  MPI_Allgather(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, srcdest, n, MPI_INT, nrnmpi_comm);
439 }
440 
441 void nrnmpi_int_allgatherv(int* s, int* r, int* n, int* dspl) {
442  MPI_Allgatherv(s, n[nrnmpi_myid], MPI_INT, r, n, dspl, MPI_INT, nrnmpi_comm);
443 }
444 
445 void nrnmpi_int_allgatherv_inplace(int* srcdest, int* n, int* dspl) {
446  MPI_Allgatherv(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, srcdest, n, dspl, MPI_INT, nrnmpi_comm);
447 }
448 
449 void nrnmpi_char_allgatherv(char* s, char* r, int* n, int* dspl) {
450  MPI_Allgatherv(s, n[nrnmpi_myid], MPI_CHAR, r, n, dspl, MPI_CHAR, nrnmpi_comm);
451 }
452 
453 void nrnmpi_long_allgatherv(int64_t* s, int64_t* r, int* n, int* dspl) {
454  MPI_Allgatherv(s, n[nrnmpi_myid], MPI_INT64_T, r, n, dspl, MPI_INT64_T, nrnmpi_comm);
455 }
456 
457 void nrnmpi_long_allgatherv_inplace(long* srcdest, int* n, int* dspl) {
458  MPI_Allgatherv(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, srcdest, n, dspl, MPI_LONG, nrnmpi_comm);
459 }
460 
461 void nrnmpi_dbl_allgatherv(double* s, double* r, int* n, int* dspl) {
462  MPI_Allgatherv(s, n[nrnmpi_myid], MPI_DOUBLE, r, n, dspl, MPI_DOUBLE, nrnmpi_comm);
463 }
464 
465 void nrnmpi_dbl_allgatherv_inplace(double* srcdest, int* n, int* dspl) {
466  MPI_Allgatherv(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, srcdest, n, dspl, MPI_DOUBLE, nrnmpi_comm);
467 }
468 
469 void nrnmpi_dbl_broadcast(double* buf, int cnt, int root) {
470  MPI_Bcast(buf, cnt, MPI_DOUBLE, root, nrnmpi_comm);
471 }
472 
473 void nrnmpi_int_broadcast(int* buf, int cnt, int root) {
474  MPI_Bcast(buf, cnt, MPI_INT, root, nrnmpi_comm);
475 }
476 
477 void nrnmpi_char_broadcast(char* buf, int cnt, int root) {
478  MPI_Bcast(buf, cnt, MPI_CHAR, root, nrnmpi_comm);
479 }
480 
481 void nrnmpi_char_broadcast_world(char** pstr, int root) {
482  int sz;
483  sz = *pstr ? (strlen(*pstr) + 1) : 0;
484  MPI_Bcast(&sz, 1, MPI_INT, root, nrnmpi_world_comm);
485  if (nrnmpi_myid_world != root) {
486  if (*pstr) {
487  free(*pstr);
488  *pstr = NULL;
489  }
490  if (sz) {
491  *pstr = static_cast<char*>(hoc_Emalloc(sz * sizeof(char)));
492  hoc_malchk();
493  }
494  }
495  if (sz) {
496  MPI_Bcast(*pstr, sz, MPI_CHAR, root, nrnmpi_world_comm);
497  }
498 }
499 
500 int nrnmpi_int_sum_reduce(int in) {
501  int result;
502  MPI_Allreduce(&in, &result, 1, MPI_INT, MPI_SUM, nrnmpi_comm);
503  return result;
504 }
505 
506 void nrnmpi_assert_opstep(int opstep, double t) {
507  /* all machines in comm should have same opstep and same t. */
508  double buf[2];
509  if (nrnmpi_numprocs < 2) {
510  return;
511  }
512  buf[0] = (double) opstep;
513  buf[1] = t;
514  MPI_Bcast(buf, 2, MPI_DOUBLE, 0, nrnmpi_comm);
515  if (opstep != (int) buf[0] || t != buf[1]) {
516  printf(
517  "%d opstep=%d %d t=%g t-troot=%g\n", nrnmpi_myid, opstep, (int) buf[0], t, t - buf[1]);
518  hoc_execerror("nrnmpi_assert_opstep failed", (char*) 0);
519  }
520 }
521 
522 double nrnmpi_dbl_allmin(double x) {
523  double result;
524  if (nrnmpi_numprocs < 2) {
525  return x;
526  }
527  MPI_Allreduce(&x, &result, 1, MPI_DOUBLE, MPI_MIN, nrnmpi_comm);
528  return result;
529 }
530 
531 static void pgvts_op(double* in, double* inout, int* len, MPI_Datatype* dptr) {
532  int i, r = 0;
533  assert(*dptr == MPI_DOUBLE);
534  assert(*len == 4);
535  if (in[0] < inout[0]) {
536  /* least time has highest priority */
537  r = 1;
538  } else if (in[0] == inout[0]) {
539  /* when times are equal then */
540  if (in[1] < inout[1]) {
541  /* NetParEvent done last */
542  r = 1;
543  } else if (in[1] == inout[1]) {
544  /* when times and ops are equal then */
545  if (in[2] < inout[2]) {
546  /* init done next to last.*/
547  r = 1;
548  } else if (in[2] == inout[2]) {
549  /* when times, ops, and inits are equal then */
550  if (in[3] < inout[3]) {
551  /* choose lowest rank */
552  r = 1;
553  }
554  }
555  }
556  }
557  if (r) {
558  for (i = 0; i < 4; ++i) {
559  inout[i] = in[i];
560  }
561  }
562 }
563 
564 int nrnmpi_pgvts_least(double* t, int* op, int* init) {
565  int i;
566  double ibuf[4], obuf[4];
567  ibuf[0] = *t;
568  ibuf[1] = (double) (*op);
569  ibuf[2] = (double) (*init);
570  ibuf[3] = (double) nrnmpi_myid;
571  for (i = 0; i < 4; ++i) {
572  obuf[i] = ibuf[i];
573  }
574  MPI_Allreduce(ibuf, obuf, 4, MPI_DOUBLE, mpi_pgvts_op, nrnmpi_comm);
575  assert(obuf[0] <= *t);
576  if (obuf[0] == *t) {
577  assert((int) obuf[1] <= *op);
578  if ((int) obuf[1] == *op) {
579  assert((int) obuf[2] <= *init);
580  if ((int) obuf[2] == *init) {
581  assert((int) obuf[3] <= nrnmpi_myid);
582  }
583  }
584  }
585  *t = obuf[0];
586  *op = (int) obuf[1];
587  *init = (int) obuf[2];
588  if (nrnmpi_myid == (int) obuf[3]) {
589  return 1;
590  }
591  return 0;
592 }
593 
594 /* following for splitcell.cpp transfer */
595 void nrnmpi_send_doubles(double* pd, int cnt, int dest, int tag) {
596  MPI_Send(pd, cnt, MPI_DOUBLE, dest, tag, nrnmpi_comm);
597 }
598 
599 void nrnmpi_recv_doubles(double* pd, int cnt, int src, int tag) {
600  MPI_Status status;
601  MPI_Recv(pd, cnt, MPI_DOUBLE, src, tag, nrnmpi_comm, &status);
602 }
603 
604 void nrnmpi_postrecv_doubles(double* pd, int cnt, int src, int tag, void** request) {
605  MPI_Irecv(pd, cnt, MPI_DOUBLE, src, tag, nrnmpi_comm, (MPI_Request*) request);
606 }
607 
608 void nrnmpi_wait(void** request) {
609  MPI_Status status;
610  MPI_Wait((MPI_Request*) request, &status);
611 }
612 
613 void nrnmpi_barrier() {
614  if (nrnmpi_numprocs < 2) {
615  return;
616  }
617  MPI_Barrier(nrnmpi_comm);
618 }
619 
620 double nrnmpi_dbl_allreduce(double x, int type) {
621  double result;
622  MPI_Op t;
623  if (nrnmpi_numprocs < 2) {
624  return x;
625  }
626  if (type == 1) {
627  t = MPI_SUM;
628  } else if (type == 2) {
629  t = MPI_MAX;
630  } else {
631  t = MPI_MIN;
632  }
633  MPI_Allreduce(&x, &result, 1, MPI_DOUBLE, t, nrnmpi_comm);
634  return result;
635 }
636 
637 extern "C" void nrnmpi_dbl_allreduce_vec(double* src, double* dest, int cnt, int type) {
638  int i;
639  MPI_Op t;
640  assert(src != dest);
641  if (nrnmpi_numprocs < 2) {
642  for (i = 0; i < cnt; ++i) {
643  dest[i] = src[i];
644  }
645  return;
646  }
647  if (type == 1) {
648  t = MPI_SUM;
649  } else if (type == 2) {
650  t = MPI_MAX;
651  } else {
652  t = MPI_MIN;
653  }
654  MPI_Allreduce(src, dest, cnt, MPI_DOUBLE, t, nrnmpi_comm);
655  return;
656 }
657 
658 void nrnmpi_longdbl_allreduce_vec(longdbl* src, longdbl* dest, int cnt, int type) {
659  int i;
660  MPI_Op t;
661  assert(src != dest);
662  if (nrnmpi_numprocs < 2) {
663  for (i = 0; i < cnt; ++i) {
664  dest[i] = src[i];
665  }
666  return;
667  }
668  if (type == 1) {
669  t = MPI_SUM;
670  } else if (type == 2) {
671  t = MPI_MAX;
672  } else {
673  t = MPI_MIN;
674  }
675  MPI_Allreduce(src, dest, cnt, MPI_LONG_DOUBLE, t, nrnmpi_comm);
676  return;
677 }
678 
679 void nrnmpi_long_allreduce_vec(long* src, long* dest, int cnt, int type) {
680  int i;
681  MPI_Op t;
682  assert(src != dest);
683  if (nrnmpi_numprocs < 2) {
684  for (i = 0; i < cnt; ++i) {
685  dest[i] = src[i];
686  }
687  return;
688  }
689  if (type == 1) {
690  t = MPI_SUM;
691  } else if (type == 2) {
692  t = MPI_MAX;
693  } else {
694  t = MPI_MIN;
695  }
696  MPI_Allreduce(src, dest, cnt, MPI_LONG, t, nrnmpi_comm);
697  return;
698 }
699 
700 void nrnmpi_dbl_allgather(double* s, double* r, int n) {
701  MPI_Allgather(s, n, MPI_DOUBLE, r, n, MPI_DOUBLE, nrnmpi_comm);
702 }
703 
704 #if BGPDMA
705 
706 static MPI_Comm bgp_comm;
707 
708 void nrnmpi_bgp_comm() {
709  if (!bgp_comm) {
710  MPI_Comm_dup(nrnmpi_comm, &bgp_comm);
711  }
712 }
713 
714 void nrnmpi_bgp_multisend(NRNMPI_Spike* spk, int n, int* hosts) {
715  int i;
716  MPI_Request r;
717  MPI_Status status;
718  for (i = 0; i < n; ++i) {
719  MPI_Isend(spk, 1, spike_type, hosts[i], 1, bgp_comm, &r);
720  MPI_Request_free(&r);
721  }
722 }
723 
725  int flag = 0;
726  MPI_Status status;
727  MPI_Iprobe(MPI_ANY_SOURCE, 1, bgp_comm, &flag, &status);
728  if (flag) {
729  MPI_Recv(spk, 1, spike_type, MPI_ANY_SOURCE, 1, bgp_comm, &status);
730  }
731  return flag;
732 }
733 
734 static int iii;
735 int nrnmpi_bgp_conserve(int nsend, int nrecv) {
736  int tcnts[2];
737  tcnts[0] = nsend - nrecv;
738  MPI_Allreduce(tcnts, tcnts + 1, 1, MPI_INT, MPI_SUM, bgp_comm);
739  return tcnts[1];
740 }
741 
742 #endif /*BGPDMA*/
743 
744 #endif /*NRNMPI*/
static void nrnmpi_int_alltoallv(int *s, int *scnt, int *sdispl, int *r, int *rcnt, int *rdispl)
static void nrnmpi_dbl_alltoallv(double *s, int *scnt, int *sdispl, double *r, int *rcnt, int *rdispl)
static void nrnmpi_dbl_allgatherv(double *s, double *r, int *n, int *dspl)
static void nrnmpi_int_allgather(int *s, int *r, int n)
static void nrnmpi_barrier()
static void nrnmpi_int_allgatherv(int *s, int *r, int *n, int *dspl)
static int nrnmpi_int_allmax(int x)
void nrnmpi_bgp_multisend(NRNMPI_Spike *, int, int *)
int nrnmpi_bgp_conserve(int nsend, int nrecv)
void nrnmpi_bgp_comm()
int nrnmpi_bgp_single_advance(NRNMPI_Spike *)
short type
Definition: cabvars.h:9
double t
Definition: cvodeobj.cpp:59
void hoc_execerror(const char *, const char *)
Definition: hoc.cpp:754
char buf[512]
Definition: init.cpp:13
#define assert(ex)
Definition: hocassrt.h:32
void nrnbbs_context_wait()
Definition: datapath.cpp:34
#define i
Definition: md1redef.h:12
void init()
Definition: init.cpp:291
int * nin_
NRNMPI_Spike * spikein_
int icapacity_
int localgid_size_
NRNMPI_Spike * spikeout_
int ovfl_capacity_
int ovfl_
unsigned char * spfixout_
int nout_
unsigned char * spfixin_
int ag_send_nspike_
unsigned char * spfixin_ovfl_
#define nrn_spikebuf_size
Definition: mpispike.h:5
int ag_send_size_
static void nrnmpi_postrecv_doubles(double *, int, int, int, void **)
Definition: multisplit.cpp:51
static void nrnmpi_wait(void **)
Definition: multisplit.cpp:53
static void nrnmpi_send_doubles(double *, int, int, int)
Definition: multisplit.cpp:52
static int nrnmpi_use
Definition: multisplit.cpp:48
#define printf
Definition: mwprefix.h:26
int const size_t const size_t n
Definition: nrngsl.h:11
return status
int nrnmpi_myid
int nrnmpi_numprocs
int nrnmpi_myid_world
MPI_Comm nrnmpi_world_comm
void * hoc_Emalloc(size_t size)
Definition: symbol.cpp:190
void hoc_malchk()
Definition: symbol.cpp:183
long double longdbl
Definition: nrnmpidec.h:10
MPI_Comm nrnmpi_comm
static void nrnmpi_dbl_broadcast(double *, int, int)
Definition: ocbbs.cpp:62
static void nrnmpi_char_broadcast(char *, int, int)
Definition: ocbbs.cpp:61
static void nrnmpi_int_broadcast(int *, int, int)
Definition: ocbbs.cpp:60
#define root
Definition: rbtqueue.cpp:53
#define cnt
Definition: spt2queue.cpp:19
#define NULL
Definition: sptree.h:16
double spiketime
Definition: nrnmpi.h:19
int gid
Definition: nrnmpi.h:18