forked from FFTW/fftw3
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmpi-bench.c
844 lines (764 loc) · 23.4 KB
/
mpi-bench.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
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
/**************************************************************************/
/* NOTE to users: this is the FFTW-MPI self-test and benchmark program.
It is probably NOT a good place to learn FFTW usage, since it has a
lot of added complexity in order to exercise and test the full API,
etcetera. We suggest reading the manual. */
/**************************************************************************/
#include <math.h>
#include <stdio.h>
#include <string.h>
#include "fftw3-mpi.h"
#include "tests/fftw-bench.h"
#if defined(BENCHFFT_SINGLE)
# define BENCH_MPI_TYPE MPI_FLOAT
#elif defined(BENCHFFT_LDOUBLE)
# define BENCH_MPI_TYPE MPI_LONG_DOUBLE
#elif defined(BENCHFFT_QUAD)
# error MPI quad-precision type is unknown
#else
# define BENCH_MPI_TYPE MPI_DOUBLE
#endif
#if SIZEOF_PTRDIFF_T == SIZEOF_INT
# define FFTW_MPI_PTRDIFF_T MPI_INT
#elif SIZEOF_PTRDIFF_T == SIZEOF_LONG
# define FFTW_MPI_PTRDIFF_T MPI_LONG
#elif SIZEOF_PTRDIFF_T == SIZEOF_LONG_LONG
# define FFTW_MPI_PTRDIFF_T MPI_LONG_LONG
#else
# error MPI type for ptrdiff_t is unknown
# define FFTW_MPI_PTRDIFF_T MPI_LONG
#endif
static const char *mkversion(void) { return FFTW(version); }
static const char *mkcc(void) { return FFTW(cc); }
static const char *mkcodelet_optim(void) { return FFTW(codelet_optim); }
static const char *mknproc(void) {
static char buf[32];
int ncpus;
MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
#ifdef HAVE_SNPRINTF
snprintf(buf, 32, "%d", ncpus);
#else
sprintf(buf, "%d", ncpus);
#endif
return buf;
}
BEGIN_BENCH_DOC
BENCH_DOC("name", "fftw3_mpi")
BENCH_DOCF("version", mkversion)
BENCH_DOCF("cc", mkcc)
BENCH_DOCF("codelet-optim", mkcodelet_optim)
BENCH_DOCF("nproc", mknproc)
END_BENCH_DOC
static int n_pes = 1, my_pe = 0;
/* global variables describing the shape of the data and its distribution */
static int rnk;
static ptrdiff_t vn, iNtot, oNtot;
static ptrdiff_t *local_ni=0, *local_starti=0;
static ptrdiff_t *local_no=0, *local_starto=0;
static ptrdiff_t *all_local_ni=0, *all_local_starti=0; /* n_pes x rnk arrays */
static ptrdiff_t *all_local_no=0, *all_local_starto=0; /* n_pes x rnk arrays */
static ptrdiff_t *istrides = 0, *ostrides = 0;
static ptrdiff_t *total_ni=0, *total_no=0;
static int *isend_cnt = 0, *isend_off = 0; /* for MPI_Scatterv */
static int *orecv_cnt = 0, *orecv_off = 0; /* for MPI_Gatherv */
static bench_real *local_in = 0, *local_out = 0;
static bench_real *all_local_in = 0, *all_local_out = 0;
static int all_local_in_alloc = 0, all_local_out_alloc = 0;
static FFTW(plan) plan_scramble_in = 0, plan_unscramble_out = 0;
static void alloc_rnk(int rnk_) {
rnk = rnk_;
bench_free(local_ni);
if (rnk == 0)
local_ni = 0;
else
local_ni = (ptrdiff_t *) bench_malloc(sizeof(ptrdiff_t) * rnk
* (8 + n_pes * 4));
local_starti = local_ni + rnk;
local_no = local_ni + 2 * rnk;
local_starto = local_ni + 3 * rnk;
istrides = local_ni + 4 * rnk;
ostrides = local_ni + 5 * rnk;
total_ni = local_ni + 6 * rnk;
total_no = local_ni + 7 * rnk;
all_local_ni = local_ni + 8 * rnk;
all_local_starti = local_ni + (8 + n_pes) * rnk;
all_local_no = local_ni + (8 + 2 * n_pes) * rnk;
all_local_starto = local_ni + (8 + 3 * n_pes) * rnk;
}
static void setup_gather_scatter(void)
{
int i, j;
ptrdiff_t off;
MPI_Gather(local_ni, rnk, FFTW_MPI_PTRDIFF_T,
all_local_ni, rnk, FFTW_MPI_PTRDIFF_T,
0, MPI_COMM_WORLD);
MPI_Bcast(all_local_ni, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
MPI_Gather(local_starti, rnk, FFTW_MPI_PTRDIFF_T,
all_local_starti, rnk, FFTW_MPI_PTRDIFF_T,
0, MPI_COMM_WORLD);
MPI_Bcast(all_local_starti, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
MPI_Gather(local_no, rnk, FFTW_MPI_PTRDIFF_T,
all_local_no, rnk, FFTW_MPI_PTRDIFF_T,
0, MPI_COMM_WORLD);
MPI_Bcast(all_local_no, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
MPI_Gather(local_starto, rnk, FFTW_MPI_PTRDIFF_T,
all_local_starto, rnk, FFTW_MPI_PTRDIFF_T,
0, MPI_COMM_WORLD);
MPI_Bcast(all_local_starto, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
off = 0;
for (i = 0; i < n_pes; ++i) {
ptrdiff_t N = vn;
for (j = 0; j < rnk; ++j)
N *= all_local_ni[i * rnk + j];
isend_cnt[i] = N;
isend_off[i] = off;
off += N;
}
iNtot = off;
all_local_in_alloc = 1;
istrides[rnk - 1] = vn;
for (j = rnk - 2; j >= 0; --j)
istrides[j] = total_ni[j + 1] * istrides[j + 1];
off = 0;
for (i = 0; i < n_pes; ++i) {
ptrdiff_t N = vn;
for (j = 0; j < rnk; ++j)
N *= all_local_no[i * rnk + j];
orecv_cnt[i] = N;
orecv_off[i] = off;
off += N;
}
oNtot = off;
all_local_out_alloc = 1;
ostrides[rnk - 1] = vn;
for (j = rnk - 2; j >= 0; --j)
ostrides[j] = total_no[j + 1] * ostrides[j + 1];
}
static void copy_block_out(const bench_real *in,
int rnk, ptrdiff_t *n, ptrdiff_t *start,
ptrdiff_t is, ptrdiff_t *os, ptrdiff_t vn,
bench_real *out)
{
ptrdiff_t i;
if (rnk == 0) {
for (i = 0; i < vn; ++i)
out[i] = in[i];
}
else if (rnk == 1) { /* this case is just an optimization */
ptrdiff_t j;
out += start[0] * os[0];
for (j = 0; j < n[0]; ++j) {
for (i = 0; i < vn; ++i)
out[i] = in[i];
in += is;
out += os[0];
}
}
else {
/* we should do n[0] for locality, but this way is simpler to code */
for (i = 0; i < n[rnk - 1]; ++i)
copy_block_out(in + i * is,
rnk - 1, n, start, is * n[rnk - 1], os, vn,
out + (start[rnk - 1] + i) * os[rnk - 1]);
}
}
static void copy_block_in(bench_real *in,
int rnk, ptrdiff_t *n, ptrdiff_t *start,
ptrdiff_t is, ptrdiff_t *os, ptrdiff_t vn,
const bench_real *out)
{
ptrdiff_t i;
if (rnk == 0) {
for (i = 0; i < vn; ++i)
in[i] = out[i];
}
else if (rnk == 1) { /* this case is just an optimization */
ptrdiff_t j;
out += start[0] * os[0];
for (j = 0; j < n[0]; ++j) {
for (i = 0; i < vn; ++i)
in[i] = out[i];
in += is;
out += os[0];
}
}
else {
/* we should do n[0] for locality, but this way is simpler to code */
for (i = 0; i < n[rnk - 1]; ++i)
copy_block_in(in + i * is,
rnk - 1, n, start, is * n[rnk - 1], os, vn,
out + (start[rnk - 1] + i) * os[rnk - 1]);
}
}
static void do_scatter_in(bench_real *in)
{
bench_real *ali;
int i;
if (all_local_in_alloc) {
bench_free(all_local_in);
all_local_in = (bench_real*) bench_malloc(iNtot*sizeof(bench_real));
all_local_in_alloc = 0;
}
ali = all_local_in;
for (i = 0; i < n_pes; ++i) {
copy_block_in(ali,
rnk, all_local_ni + i * rnk,
all_local_starti + i * rnk,
vn, istrides, vn,
in);
ali += isend_cnt[i];
}
MPI_Scatterv(all_local_in, isend_cnt, isend_off, BENCH_MPI_TYPE,
local_in, isend_cnt[my_pe], BENCH_MPI_TYPE,
0, MPI_COMM_WORLD);
}
static void do_gather_out(bench_real *out)
{
bench_real *alo;
int i;
if (all_local_out_alloc) {
bench_free(all_local_out);
all_local_out = (bench_real*) bench_malloc(oNtot*sizeof(bench_real));
all_local_out_alloc = 0;
}
MPI_Gatherv(local_out, orecv_cnt[my_pe], BENCH_MPI_TYPE,
all_local_out, orecv_cnt, orecv_off, BENCH_MPI_TYPE,
0, MPI_COMM_WORLD);
MPI_Bcast(all_local_out, oNtot, BENCH_MPI_TYPE, 0, MPI_COMM_WORLD);
alo = all_local_out;
for (i = 0; i < n_pes; ++i) {
copy_block_out(alo,
rnk, all_local_no + i * rnk,
all_local_starto + i * rnk,
vn, ostrides, vn,
out);
alo += orecv_cnt[i];
}
}
static void alloc_local(ptrdiff_t nreal, int inplace)
{
bench_free(local_in);
if (local_out != local_in) bench_free(local_out);
local_in = local_out = 0;
if (nreal > 0) {
ptrdiff_t i;
local_in = (bench_real*) bench_malloc(nreal * sizeof(bench_real));
if (inplace)
local_out = local_in;
else
local_out = (bench_real*) bench_malloc(nreal * sizeof(bench_real));
for (i = 0; i < nreal; ++i) local_in[i] = local_out[i] = 0.0;
}
}
void after_problem_rcopy_from(bench_problem *p, bench_real *ri)
{
UNUSED(p);
do_scatter_in(ri);
if (plan_scramble_in) FFTW(execute)(plan_scramble_in);
}
void after_problem_rcopy_to(bench_problem *p, bench_real *ro)
{
UNUSED(p);
if (plan_unscramble_out) FFTW(execute)(plan_unscramble_out);
do_gather_out(ro);
}
void after_problem_ccopy_from(bench_problem *p, bench_real *ri, bench_real *ii)
{
UNUSED(ii);
after_problem_rcopy_from(p, ri);
}
void after_problem_ccopy_to(bench_problem *p, bench_real *ro, bench_real *io)
{
UNUSED(io);
after_problem_rcopy_to(p, ro);
}
void after_problem_hccopy_from(bench_problem *p, bench_real *ri, bench_real *ii)
{
UNUSED(ii);
after_problem_rcopy_from(p, ri);
}
void after_problem_hccopy_to(bench_problem *p, bench_real *ro, bench_real *io)
{
UNUSED(io);
after_problem_rcopy_to(p, ro);
}
static FFTW(plan) mkplan_transpose_local(ptrdiff_t nx, ptrdiff_t ny,
ptrdiff_t vn,
bench_real *in, bench_real *out)
{
FFTW(iodim64) hdims[3];
FFTW(r2r_kind) k[3];
FFTW(plan) pln;
hdims[0].n = nx;
hdims[0].is = ny * vn;
hdims[0].os = vn;
hdims[1].n = ny;
hdims[1].is = vn;
hdims[1].os = nx * vn;
hdims[2].n = vn;
hdims[2].is = 1;
hdims[2].os = 1;
k[0] = k[1] = k[2] = FFTW_R2HC;
pln = FFTW(plan_guru64_r2r)(0, 0, 3, hdims, in, out, k, FFTW_ESTIMATE);
BENCH_ASSERT(pln != 0);
return pln;
}
static int tensor_rowmajor_transposedp(bench_tensor *t)
{
bench_iodim *d;
int i;
BENCH_ASSERT(BENCH_FINITE_RNK(t->rnk));
if (t->rnk < 2)
return 0;
d = t->dims;
if (d[0].is != d[1].is * d[1].n
|| d[0].os != d[1].is
|| d[1].os != d[0].os * d[0].n)
return 0;
if (t->rnk > 2 && d[1].is != d[2].is * d[2].n)
return 0;
for (i = 2; i + 1 < t->rnk; ++i) {
d = t->dims + i;
if (d[0].is != d[1].is * d[1].n
|| d[0].os != d[1].os * d[1].n)
return 0;
}
if (t->rnk > 2 && t->dims[t->rnk-1].is != t->dims[t->rnk-1].os)
return 0;
return 1;
}
static int tensor_contiguousp(bench_tensor *t, int s)
{
return (t->dims[t->rnk-1].is == s
&& ((tensor_rowmajorp(t) &&
t->dims[t->rnk-1].is == t->dims[t->rnk-1].os)
|| tensor_rowmajor_transposedp(t)));
}
static FFTW(plan) mkplan_complex(bench_problem *p, unsigned flags)
{
FFTW(plan) pln = 0;
int i;
ptrdiff_t ntot;
vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
if (p->sz->rnk < 1
|| p->split
|| !tensor_contiguousp(p->sz, vn)
|| tensor_rowmajor_transposedp(p->sz)
|| p->vecsz->rnk > 1
|| (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
|| p->vecsz->dims[0].os != 1)))
return 0;
alloc_rnk(p->sz->rnk);
for (i = 0; i < rnk; ++i) {
total_ni[i] = total_no[i] = p->sz->dims[i].n;
local_ni[i] = local_no[i] = total_ni[i];
local_starti[i] = local_starto[i] = 0;
}
if (rnk > 1) {
ptrdiff_t n, start, nT, startT;
ntot = FFTW(mpi_local_size_many_transposed)
(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
MPI_COMM_WORLD,
&n, &start, &nT, &startT);
if (flags & FFTW_MPI_TRANSPOSED_IN) {
local_ni[1] = nT;
local_starti[1] = startT;
}
else {
local_ni[0] = n;
local_starti[0] = start;
}
if (flags & FFTW_MPI_TRANSPOSED_OUT) {
local_no[1] = nT;
local_starto[1] = startT;
}
else {
local_no[0] = n;
local_starto[0] = start;
}
}
else if (rnk == 1) {
ntot = FFTW(mpi_local_size_many_1d)
(total_ni[0], vn, MPI_COMM_WORLD, p->sign, flags,
local_ni, local_starti, local_no, local_starto);
}
alloc_local(ntot * 2, p->in == p->out);
pln = FFTW(mpi_plan_many_dft)(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK,
FFTW_MPI_DEFAULT_BLOCK,
(FFTW(complex) *) local_in,
(FFTW(complex) *) local_out,
MPI_COMM_WORLD, p->sign, flags);
vn *= 2;
if (rnk > 1) {
ptrdiff_t nrest = 1;
for (i = 2; i < rnk; ++i) nrest *= p->sz->dims[i].n;
if (flags & FFTW_MPI_TRANSPOSED_IN)
plan_scramble_in = mkplan_transpose_local(
p->sz->dims[0].n, local_ni[1], vn * nrest,
local_in, local_in);
if (flags & FFTW_MPI_TRANSPOSED_OUT)
plan_unscramble_out = mkplan_transpose_local(
local_no[1], p->sz->dims[0].n, vn * nrest,
local_out, local_out);
}
return pln;
}
static int tensor_real_contiguousp(bench_tensor *t, int sign, int s)
{
return (t->dims[t->rnk-1].is == s
&& ((tensor_real_rowmajorp(t, sign, 1) &&
t->dims[t->rnk-1].is == t->dims[t->rnk-1].os)));
}
static FFTW(plan) mkplan_real(bench_problem *p, unsigned flags)
{
FFTW(plan) pln = 0;
int i;
ptrdiff_t ntot;
vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
if (p->sz->rnk < 2
|| p->split
|| !tensor_real_contiguousp(p->sz, p->sign, vn)
|| tensor_rowmajor_transposedp(p->sz)
|| p->vecsz->rnk > 1
|| (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
|| p->vecsz->dims[0].os != 1)))
return 0;
alloc_rnk(p->sz->rnk);
for (i = 0; i < rnk; ++i) {
total_ni[i] = total_no[i] = p->sz->dims[i].n;
local_ni[i] = local_no[i] = total_ni[i];
local_starti[i] = local_starto[i] = 0;
}
local_ni[rnk-1] = local_no[rnk-1] = total_ni[rnk-1] = total_no[rnk-1]
= p->sz->dims[rnk-1].n / 2 + 1;
{
ptrdiff_t n, start, nT, startT;
ntot = FFTW(mpi_local_size_many_transposed)
(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
MPI_COMM_WORLD,
&n, &start, &nT, &startT);
if (flags & FFTW_MPI_TRANSPOSED_IN) {
local_ni[1] = nT;
local_starti[1] = startT;
}
else {
local_ni[0] = n;
local_starti[0] = start;
}
if (flags & FFTW_MPI_TRANSPOSED_OUT) {
local_no[1] = nT;
local_starto[1] = startT;
}
else {
local_no[0] = n;
local_starto[0] = start;
}
}
alloc_local(ntot * 2, p->in == p->out);
total_ni[rnk - 1] = p->sz->dims[rnk - 1].n;
if (p->sign < 0)
pln = FFTW(mpi_plan_many_dft_r2c)(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK,
FFTW_MPI_DEFAULT_BLOCK,
local_in,
(FFTW(complex) *) local_out,
MPI_COMM_WORLD, flags);
else
pln = FFTW(mpi_plan_many_dft_c2r)(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK,
FFTW_MPI_DEFAULT_BLOCK,
(FFTW(complex) *) local_in,
local_out,
MPI_COMM_WORLD, flags);
total_ni[rnk - 1] = p->sz->dims[rnk - 1].n / 2 + 1;
vn *= 2;
{
ptrdiff_t nrest = 1;
for (i = 2; i < rnk; ++i) nrest *= total_ni[i];
if (flags & FFTW_MPI_TRANSPOSED_IN)
plan_scramble_in = mkplan_transpose_local(
total_ni[0], local_ni[1], vn * nrest,
local_in, local_in);
if (flags & FFTW_MPI_TRANSPOSED_OUT)
plan_unscramble_out = mkplan_transpose_local(
local_no[1], total_ni[0], vn * nrest,
local_out, local_out);
}
return pln;
}
static FFTW(plan) mkplan_transpose(bench_problem *p, unsigned flags)
{
ptrdiff_t ntot, nx, ny;
int ix=0, iy=1, i;
const bench_iodim *d = p->vecsz->dims;
FFTW(plan) pln;
if (p->vecsz->rnk == 3) {
for (i = 0; i < 3; ++i)
if (d[i].is == 1 && d[i].os == 1) {
vn = d[i].n;
ix = (i + 1) % 3;
iy = (i + 2) % 3;
break;
}
if (i == 3) return 0;
}
else {
vn = 1;
ix = 0;
iy = 1;
}
if (d[ix].is == d[iy].n * vn && d[ix].os == vn
&& d[iy].os == d[ix].n * vn && d[iy].is == vn) {
nx = d[ix].n;
ny = d[iy].n;
}
else if (d[iy].is == d[ix].n * vn && d[iy].os == vn
&& d[ix].os == d[iy].n * vn && d[ix].is == vn) {
nx = d[iy].n;
ny = d[ix].n;
}
else
return 0;
alloc_rnk(2);
ntot = vn * FFTW(mpi_local_size_2d_transposed)(nx, ny, MPI_COMM_WORLD,
&local_ni[0],
&local_starti[0],
&local_no[0],
&local_starto[0]);
local_ni[1] = ny;
local_starti[1] = 0;
local_no[1] = nx;
local_starto[1] = 0;
total_ni[0] = nx; total_ni[1] = ny;
total_no[1] = nx; total_no[0] = ny;
alloc_local(ntot, p->in == p->out);
pln = FFTW(mpi_plan_many_transpose)(nx, ny, vn,
FFTW_MPI_DEFAULT_BLOCK,
FFTW_MPI_DEFAULT_BLOCK,
local_in, local_out,
MPI_COMM_WORLD, flags);
if (flags & FFTW_MPI_TRANSPOSED_IN)
plan_scramble_in = mkplan_transpose_local(local_ni[0], ny, vn,
local_in, local_in);
if (flags & FFTW_MPI_TRANSPOSED_OUT)
plan_unscramble_out = mkplan_transpose_local
(nx, local_no[0], vn, local_out, local_out);
#if 0
if (pln && vn == 1) {
int i, j;
bench_real *ri = (bench_real *) p->in;
bench_real *ro = (bench_real *) p->out;
if (!ri || !ro) return pln;
setup_gather_scatter();
for (i = 0; i < nx * ny; ++i)
ri[i] = i;
after_problem_rcopy_from(p, ri);
FFTW(execute)(pln);
after_problem_rcopy_to(p, ro);
if (my_pe == 0) {
for (i = 0; i < nx; ++i) {
for (j = 0; j < ny; ++j)
printf(" %3g", ro[j * nx + i]);
printf("\n");
}
}
}
#endif
return pln;
}
static FFTW(plan) mkplan_r2r(bench_problem *p, unsigned flags)
{
FFTW(plan) pln = 0;
int i;
ptrdiff_t ntot;
FFTW(r2r_kind) *k;
if ((p->sz->rnk == 0 || (p->sz->rnk == 1 && p->sz->dims[0].n == 1))
&& p->vecsz->rnk >= 2 && p->vecsz->rnk <= 3)
return mkplan_transpose(p, flags);
vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
if (p->sz->rnk < 1
|| p->split
|| !tensor_contiguousp(p->sz, vn)
|| tensor_rowmajor_transposedp(p->sz)
|| p->vecsz->rnk > 1
|| (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
|| p->vecsz->dims[0].os != 1)))
return 0;
alloc_rnk(p->sz->rnk);
for (i = 0; i < rnk; ++i) {
total_ni[i] = total_no[i] = p->sz->dims[i].n;
local_ni[i] = local_no[i] = total_ni[i];
local_starti[i] = local_starto[i] = 0;
}
if (rnk > 1) {
ptrdiff_t n, start, nT, startT;
ntot = FFTW(mpi_local_size_many_transposed)
(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
MPI_COMM_WORLD,
&n, &start, &nT, &startT);
if (flags & FFTW_MPI_TRANSPOSED_IN) {
local_ni[1] = nT;
local_starti[1] = startT;
}
else {
local_ni[0] = n;
local_starti[0] = start;
}
if (flags & FFTW_MPI_TRANSPOSED_OUT) {
local_no[1] = nT;
local_starto[1] = startT;
}
else {
local_no[0] = n;
local_starto[0] = start;
}
}
else if (rnk == 1) {
ntot = FFTW(mpi_local_size_many_1d)
(total_ni[0], vn, MPI_COMM_WORLD, p->sign, flags,
local_ni, local_starti, local_no, local_starto);
}
alloc_local(ntot, p->in == p->out);
k = (FFTW(r2r_kind) *) bench_malloc(sizeof(FFTW(r2r_kind)) * p->sz->rnk);
for (i = 0; i < p->sz->rnk; ++i)
switch (p->k[i]) {
case R2R_R2HC: k[i] = FFTW_R2HC; break;
case R2R_HC2R: k[i] = FFTW_HC2R; break;
case R2R_DHT: k[i] = FFTW_DHT; break;
case R2R_REDFT00: k[i] = FFTW_REDFT00; break;
case R2R_REDFT01: k[i] = FFTW_REDFT01; break;
case R2R_REDFT10: k[i] = FFTW_REDFT10; break;
case R2R_REDFT11: k[i] = FFTW_REDFT11; break;
case R2R_RODFT00: k[i] = FFTW_RODFT00; break;
case R2R_RODFT01: k[i] = FFTW_RODFT01; break;
case R2R_RODFT10: k[i] = FFTW_RODFT10; break;
case R2R_RODFT11: k[i] = FFTW_RODFT11; break;
default: BENCH_ASSERT(0);
}
pln = FFTW(mpi_plan_many_r2r)(p->sz->rnk, total_ni, vn,
FFTW_MPI_DEFAULT_BLOCK,
FFTW_MPI_DEFAULT_BLOCK,
local_in, local_out,
MPI_COMM_WORLD, k, flags);
bench_free(k);
if (rnk > 1) {
ptrdiff_t nrest = 1;
for (i = 2; i < rnk; ++i) nrest *= p->sz->dims[i].n;
if (flags & FFTW_MPI_TRANSPOSED_IN)
plan_scramble_in = mkplan_transpose_local(
p->sz->dims[0].n, local_ni[1], vn * nrest,
local_in, local_in);
if (flags & FFTW_MPI_TRANSPOSED_OUT)
plan_unscramble_out = mkplan_transpose_local(
local_no[1], p->sz->dims[0].n, vn * nrest,
local_out, local_out);
}
return pln;
}
FFTW(plan) mkplan(bench_problem *p, unsigned flags)
{
FFTW(plan) pln = 0;
FFTW(destroy_plan)(plan_scramble_in); plan_scramble_in = 0;
FFTW(destroy_plan)(plan_unscramble_out); plan_unscramble_out = 0;
if (p->scrambled_in) {
if (p->sz->rnk == 1 && p->sz->dims[0].n != 1)
flags |= FFTW_MPI_SCRAMBLED_IN;
else
flags |= FFTW_MPI_TRANSPOSED_IN;
}
if (p->scrambled_out) {
if (p->sz->rnk == 1 && p->sz->dims[0].n != 1)
flags |= FFTW_MPI_SCRAMBLED_OUT;
else
flags |= FFTW_MPI_TRANSPOSED_OUT;
}
switch (p->kind) {
case PROBLEM_COMPLEX:
pln =mkplan_complex(p, flags);
break;
case PROBLEM_REAL:
pln = mkplan_real(p, flags);
break;
case PROBLEM_R2R:
pln = mkplan_r2r(p, flags);
break;
default: BENCH_ASSERT(0);
}
if (pln) setup_gather_scatter();
return pln;
}
void main_init(int *argc, char ***argv)
{
#ifdef HAVE_SMP
# if MPI_VERSION >= 2 /* for MPI_Init_thread */
int provided;
MPI_Init_thread(argc, argv, MPI_THREAD_FUNNELED, &provided);
threads_ok = provided >= MPI_THREAD_FUNNELED;
# else
MPI_Init(argc, argv);
threads_ok = 0;
# endif
#else
MPI_Init(argc, argv);
#endif
MPI_Comm_rank(MPI_COMM_WORLD, &my_pe);
MPI_Comm_size(MPI_COMM_WORLD, &n_pes);
if (my_pe != 0) verbose = -999;
no_speed_allocation = 1; /* so we can benchmark transforms > memory */
always_pad_real = 1; /* out-of-place real transforms are padded */
isend_cnt = (int *) bench_malloc(sizeof(int) * n_pes);
isend_off = (int *) bench_malloc(sizeof(int) * n_pes);
orecv_cnt = (int *) bench_malloc(sizeof(int) * n_pes);
orecv_off = (int *) bench_malloc(sizeof(int) * n_pes);
/* init_threads must be called before any other FFTW function,
including mpi_init, because it has to register the threads hooks
before the planner is initalized */
#ifdef HAVE_SMP
if (threads_ok) { BENCH_ASSERT(FFTW(init_threads)()); }
#endif
FFTW(mpi_init)();
}
void initial_cleanup(void)
{
alloc_rnk(0);
alloc_local(0, 0);
bench_free(all_local_in); all_local_in = 0;
bench_free(all_local_out); all_local_out = 0;
bench_free(isend_off); isend_off = 0;
bench_free(isend_cnt); isend_cnt = 0;
bench_free(orecv_off); orecv_off = 0;
bench_free(orecv_cnt); orecv_cnt = 0;
FFTW(destroy_plan)(plan_scramble_in); plan_scramble_in = 0;
FFTW(destroy_plan)(plan_unscramble_out); plan_unscramble_out = 0;
}
void final_cleanup(void)
{
MPI_Finalize();
}
void bench_exit(int status)
{
MPI_Abort(MPI_COMM_WORLD, status);
}
double bench_cost_postprocess(double cost)
{
double cost_max;
MPI_Allreduce(&cost, &cost_max, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
return cost_max;
}
int import_wisdom(FILE *f)
{
int success = 1, sall;
if (my_pe == 0) success = FFTW(import_wisdom_from_file)(f);
FFTW(mpi_broadcast_wisdom)(MPI_COMM_WORLD);
MPI_Allreduce(&success, &sall, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
return sall;
}
void export_wisdom(FILE *f)
{
FFTW(mpi_gather_wisdom)(MPI_COMM_WORLD);
if (my_pe == 0) FFTW(export_wisdom_to_file)(f);
}