-
Notifications
You must be signed in to change notification settings - Fork 12
/
collcomm.c
685 lines (584 loc) · 21 KB
/
collcomm.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
/***********************************************************************/
/* */
/* The Caml/MPI interface */
/* */
/* Xavier Leroy, projet Cristal, INRIA Rocquencourt */
/* */
/* Copyright 1998 Institut National de Recherche en Informatique et */
/* en Automatique. All rights reserved. This file is distributed */
/* under the terms of the GNU Library General Public License, with */
/* the special exception on linking described in file LICENSE. */
/* */
/***********************************************************************/
/* $Id$ */
/* Group communication */
#include <mpi.h>
#include <caml/mlvalues.h>
#include <caml/memory.h>
#include <caml/alloc.h>
#include <caml/bigarray.h>
#include "camlmpi.h"
/* Barrier synchronization */
value caml_mpi_barrier(value comm)
{
MPI_Barrier(Comm_val(comm));
return Val_unit;
}
/* Broadcast */
value caml_mpi_broadcast(value buffer, value root, value comm)
{
MPI_Bcast(Bp_val(buffer), caml_string_length(buffer), MPI_BYTE,
Int_val(root), Comm_val(comm));
return Val_unit;
}
value caml_mpi_broadcast_int(value data, value root, value comm)
{
long n = Long_val(data);
MPI_Bcast(&n, 1, MPI_LONG, Int_val(root), Comm_val(comm));
return Val_long(n);
}
value caml_mpi_broadcast_float(value data, value root, value comm)
{
double d = Double_val(data);
MPI_Bcast(&d, 1, MPI_DOUBLE, Int_val(root), Comm_val(comm));
return caml_copy_double(d);
}
value caml_mpi_broadcast_intarray(value data, value root, value comm)
{
MPI_Bcast(Longptr_val(data), Wosize_val(data), MPI_LONG,
Int_val(root), Comm_val(comm));
return Val_unit;
}
value caml_mpi_broadcast_floatarray(value data, value root, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
MPI_Bcast(d, len, MPI_DOUBLE, Int_val(root), Comm_val(comm));
caml_mpi_commit_floatarray(d, data, len);
return Val_unit;
}
value caml_mpi_broadcast_bigarray(value data, value root, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[d->flags & CAML_BA_KIND_MASK];
MPI_Bcast(d->data, dlen, dt, Int_val(root), Comm_val(comm));
return Val_unit;
}
/* Scatter */
static void caml_mpi_counts_displs(value lengths,
/* out */ int ** counts,
/* out */ int ** displs)
{
int size, disp, i;
size = Wosize_val(lengths);
if (size > 0) {
*counts = caml_stat_alloc(size * sizeof(int));
*displs = caml_stat_alloc(size * sizeof(int));
for (i = 0, disp = 0; i < size; i++) {
(*counts)[i] = Int_val(Field(lengths, i));
(*displs)[i] = disp;
disp += (*counts)[i];
}
} else {
*counts = NULL;
*displs = NULL;
}
}
value caml_mpi_scatter(value sendbuf, value sendlengths,
value recvbuf,
value root, value comm)
{
int * sendcounts, * displs;
caml_mpi_counts_displs(sendlengths, &sendcounts, &displs);
MPI_Scatterv(String_val(sendbuf), sendcounts, displs, MPI_BYTE,
Bp_val(recvbuf), caml_string_length(recvbuf), MPI_BYTE,
Int_val(root), Comm_val(comm));
if (sendcounts != NULL) {
caml_stat_free(sendcounts);
caml_stat_free(displs);
}
return Val_unit;
}
value caml_mpi_scatter_int(value data, value root, value comm)
{
value n;
MPI_Scatter(Longptr_val(data), 1, MPI_LONG,
&n, 1, MPI_LONG,
Int_val(root), Comm_val(comm));
return n;
}
value caml_mpi_scatter_float(value data, value root, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
double * src = caml_mpi_input_floatarray(data, len);
double dst;
MPI_Scatter(src, 1, MPI_DOUBLE, &dst, 1, MPI_DOUBLE,
Int_val(root), Comm_val(comm));
caml_mpi_free_floatarray(src);
return caml_copy_double(dst);
}
CAMLprim value caml_mpi_scatter_from_bigarray(value data, value root,
value comm)
{
CAMLparam3(data, root, comm);
struct caml_ba_array* d = Caml_ba_array_val(data);
intnat kind = d->flags & CAML_BA_KIND_MASK;
MPI_Comm c = Comm_val(comm);
int rank, csize;
MPI_Datatype dt = caml_mpi_ba_mpi_type[kind];
any_ba_value(dst);
MPI_Comm_rank(c, &rank);
if (rank == Int_val(root)) {
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(d) != csize)
caml_mpi_raise_error("Mpi.scatter_from_bigarray: array size mismatch");
}
MPI_Scatter(d->data, 1, dt, &dst, 1, dt, Int_val(root), c);
CAMLreturn(caml_mpi_ba_value(dst, kind));
}
value caml_mpi_scatter_intarray(value source, value dest,
value root, value comm)
{
mlsize_t len = Wosize_val(dest);
MPI_Scatter(Longptr_val(source), len, MPI_LONG,
Longptr_val(dest), len, MPI_LONG,
Int_val(root), Comm_val(comm));
return Val_unit;
}
value caml_mpi_scatter_floatarray(value source, value dest,
value root, value comm)
{
mlsize_t srclen = Wosize_val(source) / Double_wosize;
mlsize_t len = Wosize_val(dest) / Double_wosize;
double * src = caml_mpi_input_floatarray_at_node(source, srclen, root, comm);
double * dst = caml_mpi_output_floatarray(dest, len);
MPI_Scatter(src, len, MPI_DOUBLE, dst, len, MPI_DOUBLE,
Int_val(root), Comm_val(comm));
caml_mpi_free_floatarray(src);
caml_mpi_commit_floatarray(dst, dest, len);
return Val_unit;
}
value caml_mpi_scatter_bigarray(value source, value dest,
value root, value comm)
{
struct caml_ba_array* s = Caml_ba_array_val(source);
struct caml_ba_array* d = Caml_ba_array_val(dest);
MPI_Comm c = Comm_val(comm);
int rank, csize;
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[d->flags & CAML_BA_KIND_MASK];
MPI_Comm_rank(c, &rank);
if (rank == Int_val(root)) {
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(s) != dlen * csize)
caml_mpi_raise_error("Mpi.scatter_bigarray: array size mismatch");
}
MPI_Scatter(s->data, dlen, dt, d->data, dlen, dt, Int_val(root), c);
return Val_unit;
}
/* Gather */
value caml_mpi_gather(value sendbuf,
value recvbuf, value recvlengths,
value root, value comm)
{
int * recvcounts, * displs;
caml_mpi_counts_displs(recvlengths, &recvcounts, &displs);
MPI_Gatherv(String_val(sendbuf), caml_string_length(sendbuf), MPI_BYTE,
Bp_val(recvbuf), recvcounts, displs, MPI_BYTE,
Int_val(root), Comm_val(comm));
if (recvcounts != NULL) {
caml_stat_free(recvcounts);
caml_stat_free(displs);
}
return Val_unit;
}
value caml_mpi_gather_int(value data, value result, value root, value comm)
{
MPI_Gather(&data, 1, MPI_LONG,
Longptr_val(result), 1, MPI_LONG,
Int_val(root), Comm_val(comm));
return Val_unit;
}
value caml_mpi_gather_intarray(value data, value result,
value root, value comm)
{
mlsize_t len = Wosize_val(data);
MPI_Gather(Longptr_val(data), len, MPI_LONG,
Longptr_val(result), len, MPI_LONG,
Int_val(root), Comm_val(comm));
return Val_unit;
}
value caml_mpi_gather_float(value data, value result, value root, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
mlsize_t reslen = Wosize_val(result) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res =
caml_mpi_output_floatarray_at_node(result, reslen, root, comm);
MPI_Gather(d, len, MPI_DOUBLE, res, len, MPI_DOUBLE,
Int_val(root), Comm_val(comm));
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, reslen);
return Val_unit;
}
CAMLprim value caml_mpi_gather_to_bigarray(value data, value result,
value root, value comm)
{
CAMLparam4(data, result, root, comm);
struct caml_ba_array* r = Caml_ba_array_val(result);
intnat kind = r->flags & CAML_BA_KIND_MASK;
MPI_Comm c = Comm_val(comm);
int rank, csize;
MPI_Datatype dt = caml_mpi_ba_mpi_type[kind];
any_ba_value(d);
MPI_Comm_rank(c, &rank);
if (rank == Int_val(root)) {
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(r) != csize)
caml_mpi_raise_error("Mpi.gather_to_bigarray: array size mismatch");
}
caml_mpi_ba_element(data, kind, d);
MPI_Gather(d, 1, dt, r->data, 1, dt, Int_val(root), c);
CAMLreturn(Val_unit);
}
value caml_mpi_gather_bigarray(value data, value result,
value root, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
MPI_Comm c = Comm_val(comm);
int rank, csize;
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[r->flags & CAML_BA_KIND_MASK];
MPI_Comm_rank(c, &rank);
if (rank == Int_val(root)) {
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(r) != dlen * csize)
caml_mpi_raise_error("Mpi.gather_bigarray: array size mismatch");
}
MPI_Gather(d->data, dlen, dt, r->data, dlen, dt, Int_val(root), c);
return Val_unit;
}
/* Gather to all */
value caml_mpi_allgather(value sendbuf,
value recvbuf, value recvlengths,
value comm)
{
int * recvcounts, * displs;
caml_mpi_counts_displs(recvlengths, &recvcounts, &displs);
MPI_Allgatherv(String_val(sendbuf), caml_string_length(sendbuf), MPI_BYTE,
Bp_val(recvbuf), recvcounts, displs, MPI_BYTE,
Comm_val(comm));
caml_stat_free(recvcounts);
caml_stat_free(displs);
return Val_unit;
}
value caml_mpi_allgather_int(value data, value result, value comm)
{
MPI_Allgather(&data, 1, MPI_LONG,
Longptr_val(result), 1, MPI_LONG,
Comm_val(comm));
return Val_unit;
}
value caml_mpi_allgather_intarray(value data, value result, value comm)
{
mlsize_t len = Wosize_val(data);
MPI_Allgather(Longptr_val(data), len, MPI_LONG,
Longptr_val(result), len, MPI_LONG,
Comm_val(comm));
return Val_unit;
}
value caml_mpi_allgather_float(value data, value result, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
mlsize_t reslen = Wosize_val(result) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res = caml_mpi_output_floatarray(result, reslen);
MPI_Allgather(d, len, MPI_DOUBLE, res, len, MPI_DOUBLE,
Comm_val(comm));
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, reslen);
return Val_unit;
}
CAMLprim value caml_mpi_allgather_to_bigarray(value data, value result,
value comm)
{
CAMLparam3(data, result, comm);
struct caml_ba_array* r = Caml_ba_array_val(result);
intnat kind = r->flags & CAML_BA_KIND_MASK;
MPI_Comm c = Comm_val(comm);
int csize;
MPI_Datatype dt = caml_mpi_ba_mpi_type[kind];
any_ba_value(d);
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(r) != csize)
caml_mpi_raise_error("Mpi.allgather_to_bigarray: array size mismatch");
caml_mpi_ba_element(data, kind, d);
MPI_Allgather(d, 1, dt, r->data, 1, dt, c);
CAMLreturn(Val_unit);
}
value caml_mpi_allgather_bigarray(value data, value result, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
MPI_Comm c = Comm_val(comm);
int csize;
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[r->flags & CAML_BA_KIND_MASK];
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(r) != dlen * csize)
caml_mpi_raise_error("Mpi.allgather_bigarray: array size mismatch");
MPI_Allgather(d->data, dlen, dt, r->data, dlen, dt, c);
return Val_unit;
}
/* All to all */
value caml_mpi_alltoall(value sendbuf, value sendlengths,
value recvbuf, value recvlengths,
value comm)
{
int * recvcounts, * recvdispls;
int * sendcounts, * senddispls;
caml_mpi_counts_displs(sendlengths, &sendcounts, &senddispls);
caml_mpi_counts_displs(recvlengths, &recvcounts, &recvdispls);
MPI_Alltoallv(String_val(sendbuf), sendcounts, senddispls, MPI_BYTE,
Bp_val(recvbuf), recvcounts, recvdispls, MPI_BYTE,
Comm_val(comm));
caml_stat_free(recvcounts);
caml_stat_free(recvdispls);
caml_stat_free(sendcounts);
caml_stat_free(senddispls);
return Val_unit;
}
value caml_mpi_alltoall_intarray(value data, value result, value comm)
{
mlsize_t len = Wosize_val(data);
MPI_Comm c = Comm_val(comm);
int csize, count;
void* sendbuf = Longptr_val(data);
void* recvbuf = Longptr_val(result);
MPI_Comm_size(c, &csize);
count = len / csize;
if (len % csize != 0)
caml_mpi_raise_error("Mpi.alltoall_intarray: incorrect array size");
if (sendbuf == recvbuf) sendbuf = MPI_IN_PLACE;
MPI_Alltoall(sendbuf, count, MPI_LONG, recvbuf, count, MPI_LONG, c);
return Val_unit;
}
value caml_mpi_alltoall_floatarray(value data, value result, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
mlsize_t reslen = Wosize_val(result) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res = caml_mpi_output_floatarray(result, reslen);
MPI_Comm c = Comm_val(comm);
int csize, count;
MPI_Comm_size(c, &csize);
count = len / csize;
if (len % csize != 0)
caml_mpi_raise_error("Mpi.alltoall_floatarray: incorrect array size");
MPI_Alltoall(d, count, MPI_DOUBLE, res, count, MPI_DOUBLE, c);
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, reslen);
return Val_unit;
}
value caml_mpi_alltoall_bigarray(value data, value result, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
MPI_Comm c = Comm_val(comm);
int csize, count;
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[r->flags & CAML_BA_KIND_MASK];
MPI_Comm_size(c, &csize);
if (caml_ba_num_elts(r) != dlen)
caml_mpi_raise_error("Mpi.alltoall_bigarray: array size mismatch");
count = dlen / csize;
if (dlen % csize != 0)
caml_mpi_raise_error("Mpi.alltoall_bigarray: incorrect array size");
MPI_Alltoall(d->data, count, dt, r->data, count, dt, c);
return Val_unit;
}
/* Reduce */
static MPI_Op reduce_op[] =
{ MPI_MAX, MPI_MIN, MPI_SUM, MPI_PROD, MPI_BAND, MPI_BOR, MPI_BXOR,
// deprecated Int_* ops:
MPI_MAX, MPI_MIN, MPI_SUM, MPI_PROD, MPI_BAND, MPI_BOR, MPI_BXOR,
// deprecated Float_* ops:
MPI_MAX, MPI_MIN, MPI_SUM, MPI_PROD };
value caml_mpi_reduce_int(value data, value op, value root, value comm)
{
long d = Long_val(data);
long r = 0;
MPI_Reduce(&d, &r, 1, MPI_LONG,
reduce_op[Int_val(op)], Int_val(root), Comm_val(comm));
return Val_long(r);
}
value caml_mpi_reduce_intarray(value data, value result, value op,
value root, value comm)
{
mlsize_t len = Wosize_val(data);
int myrank;
/* Decode data at all nodes in place */
caml_mpi_decode_intarray(data, len);
/* Do the reduce */
MPI_Reduce(Longptr_val(data), Longptr_val(result), len, MPI_LONG,
reduce_op[Int_val(op)], Int_val(root), Comm_val(comm));
/* Re-encode data at all nodes in place */
caml_mpi_encode_intarray(data, len);
/* At root node, also encode result */
MPI_Comm_rank(Comm_val(comm), &myrank);
if (myrank == Int_val(root)) caml_mpi_encode_intarray(result, len);
return Val_unit;
}
value caml_mpi_reduce_float(value data, value op, value root, value comm)
{
double d = Double_val(data);
double r = 0.0;
MPI_Reduce(&d, &r, 1, MPI_DOUBLE,
reduce_op[Int_val(op)], Int_val(root), Comm_val(comm));
return caml_copy_double(r);
}
value caml_mpi_reduce_floatarray(value data, value result, value op,
value root, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res = caml_mpi_output_floatarray(result, len);
MPI_Reduce(d, res, len, MPI_DOUBLE,
reduce_op[Int_val(op)], Int_val(root), Comm_val(comm));
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, len);
return Val_unit;
}
value caml_mpi_reduce_bigarray(value data, value result, value op,
value root, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
MPI_Comm c = Comm_val(comm);
int rank;
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[d->flags & CAML_BA_KIND_MASK];
void* sendbuf = d->data;
MPI_Comm_rank(c, &rank);
if (rank == root) {
if (dlen != caml_ba_num_elts(r))
caml_mpi_raise_error("Mpi.reduce_bigarray: array size mismatch");
if (d->data == r->data) sendbuf = MPI_IN_PLACE;
}
MPI_Reduce(sendbuf, r->data, dlen, dt,
reduce_op[Int_val(op)], Int_val(root), c);
return Val_unit;
}
/* Allreduce */
value caml_mpi_allreduce_int(value data, value op, value comm)
{
long d = Long_val(data);
long r;
MPI_Allreduce(&d, &r, 1, MPI_LONG,
reduce_op[Int_val(op)], Comm_val(comm));
return Val_long(r);
}
value caml_mpi_allreduce_intarray(value data, value result, value op,
value comm)
{
mlsize_t len = Wosize_val(data);
/* Decode data at all nodes in place */
caml_mpi_decode_intarray(data, len);
/* Do the reduce */
MPI_Allreduce(Longptr_val(data), Longptr_val(result), len, MPI_LONG,
reduce_op[Int_val(op)], Comm_val(comm));
/* Re-encode data at all nodes in place */
caml_mpi_encode_intarray(data, len);
/* Re-encode result at all nodes in place */
caml_mpi_encode_intarray(result, len);
return Val_unit;
}
value caml_mpi_allreduce_float(value data, value op, value comm)
{
double d = Double_val(data);
double r;
MPI_Allreduce(&d, &r, 1, MPI_DOUBLE,
reduce_op[Int_val(op)], Comm_val(comm));
return caml_copy_double(r);
}
value caml_mpi_allreduce_floatarray(value data, value result, value op,
value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res = caml_mpi_output_floatarray(result, len);
MPI_Allreduce(d, res, len, MPI_DOUBLE,
reduce_op[Int_val(op)], Comm_val(comm));
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, len);
return Val_unit;
}
value caml_mpi_allreduce_bigarray(value data, value result, value op,
value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[d->flags & CAML_BA_KIND_MASK];
void* sendbuf = (d->data == r->data) ? MPI_IN_PLACE : d->data;
if (caml_ba_num_elts(r) != dlen)
caml_mpi_raise_error("Mpi.allreduce_bigarray: array size mismatch");
MPI_Allreduce(sendbuf, r->data, dlen, dt,
reduce_op[Int_val(op)], Comm_val(comm));
return Val_unit;
}
/* Scan */
value caml_mpi_scan_int(value data, value op, value comm)
{
long d = Long_val(data);
long r;
MPI_Scan(&d, &r, 1, MPI_LONG, reduce_op[Int_val(op)], Comm_val(comm));
return Val_long(r);
}
value caml_mpi_scan_intarray(value data, value result, value op, value comm)
{
mlsize_t len = Wosize_val(data);
/* Decode data at all nodes in place */
caml_mpi_decode_intarray(data, len);
/* Do the scan */
MPI_Scan(Longptr_val(data), Longptr_val(result), len, MPI_LONG,
reduce_op[Int_val(op)], Comm_val(comm));
/* Re-encode data at all nodes in place */
caml_mpi_encode_intarray(data, len);
/* Encode result */
caml_mpi_encode_intarray(result, len);
return Val_unit;
}
value caml_mpi_scan_float(value data, value op, value comm)
{
double d = Double_val(data), r;
MPI_Scan(&d, &r, 1, MPI_DOUBLE,
reduce_op[Int_val(op)], Comm_val(comm));
return caml_copy_double(r);
}
value caml_mpi_scan_floatarray(value data, value result, value op, value comm)
{
mlsize_t len = Wosize_val(data) / Double_wosize;
double * d = caml_mpi_input_floatarray(data, len);
double * res = caml_mpi_output_floatarray(result, len);
MPI_Scan(d, res, len, MPI_DOUBLE,
reduce_op[Int_val(op)], Comm_val(comm));
caml_mpi_free_floatarray(d);
caml_mpi_commit_floatarray(res, result, len);
return Val_unit;
}
value caml_mpi_scan_bigarray(value data, value result, value op, value comm)
{
struct caml_ba_array* d = Caml_ba_array_val(data);
struct caml_ba_array* r = Caml_ba_array_val(result);
mlsize_t dlen = caml_ba_num_elts(d);
MPI_Datatype dt = caml_mpi_ba_mpi_type[d->flags & CAML_BA_KIND_MASK];
void* sendbuf = (d->data == r->data) ? MPI_IN_PLACE : d->data;
if (caml_ba_num_elts(r) != dlen)
caml_mpi_raise_error("Mpi.scan_bigarray: array size mismatch");
MPI_Scan(sendbuf, r->data, dlen, dt,
reduce_op[Int_val(op)], Comm_val(comm));
return Val_unit;
}