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