annotate src/fftw-3.3.5/mpi/mpi-bench.c @ 43:5ea0608b923f

Current zlib source
author Chris Cannam
date Tue, 18 Oct 2016 14:33:52 +0100
parents 2cd0e3b3e1fd
children
rev   line source
Chris@42 1 /**************************************************************************/
Chris@42 2 /* NOTE to users: this is the FFTW-MPI self-test and benchmark program.
Chris@42 3 It is probably NOT a good place to learn FFTW usage, since it has a
Chris@42 4 lot of added complexity in order to exercise and test the full API,
Chris@42 5 etcetera. We suggest reading the manual. */
Chris@42 6 /**************************************************************************/
Chris@42 7
Chris@42 8 #include <math.h>
Chris@42 9 #include <stdio.h>
Chris@42 10 #include <string.h>
Chris@42 11 #include "fftw3-mpi.h"
Chris@42 12 #include "fftw-bench.h"
Chris@42 13
Chris@42 14 #if defined(BENCHFFT_SINGLE)
Chris@42 15 # define BENCH_MPI_TYPE MPI_FLOAT
Chris@42 16 #elif defined(BENCHFFT_LDOUBLE)
Chris@42 17 # define BENCH_MPI_TYPE MPI_LONG_DOUBLE
Chris@42 18 #elif defined(BENCHFFT_QUAD)
Chris@42 19 # error MPI quad-precision type is unknown
Chris@42 20 #else
Chris@42 21 # define BENCH_MPI_TYPE MPI_DOUBLE
Chris@42 22 #endif
Chris@42 23
Chris@42 24 #if SIZEOF_PTRDIFF_T == SIZEOF_INT
Chris@42 25 # define FFTW_MPI_PTRDIFF_T MPI_INT
Chris@42 26 #elif SIZEOF_PTRDIFF_T == SIZEOF_LONG
Chris@42 27 # define FFTW_MPI_PTRDIFF_T MPI_LONG
Chris@42 28 #elif SIZEOF_PTRDIFF_T == SIZEOF_LONG_LONG
Chris@42 29 # define FFTW_MPI_PTRDIFF_T MPI_LONG_LONG
Chris@42 30 #else
Chris@42 31 # error MPI type for ptrdiff_t is unknown
Chris@42 32 # define FFTW_MPI_PTRDIFF_T MPI_LONG
Chris@42 33 #endif
Chris@42 34
Chris@42 35 static const char *mkversion(void) { return FFTW(version); }
Chris@42 36 static const char *mkcc(void) { return FFTW(cc); }
Chris@42 37 static const char *mkcodelet_optim(void) { return FFTW(codelet_optim); }
Chris@42 38 static const char *mknproc(void) {
Chris@42 39 static char buf[32];
Chris@42 40 int ncpus;
Chris@42 41 MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
Chris@42 42 #ifdef HAVE_SNPRINTF
Chris@42 43 snprintf(buf, 32, "%d", ncpus);
Chris@42 44 #else
Chris@42 45 sprintf(buf, "%d", ncpus);
Chris@42 46 #endif
Chris@42 47 return buf;
Chris@42 48 }
Chris@42 49
Chris@42 50 BEGIN_BENCH_DOC
Chris@42 51 BENCH_DOC("name", "fftw3_mpi")
Chris@42 52 BENCH_DOCF("version", mkversion)
Chris@42 53 BENCH_DOCF("cc", mkcc)
Chris@42 54 BENCH_DOCF("codelet-optim", mkcodelet_optim)
Chris@42 55 BENCH_DOCF("nproc", mknproc)
Chris@42 56 END_BENCH_DOC
Chris@42 57
Chris@42 58 static int n_pes = 1, my_pe = 0;
Chris@42 59
Chris@42 60 /* global variables describing the shape of the data and its distribution */
Chris@42 61 static int rnk;
Chris@42 62 static ptrdiff_t vn, iNtot, oNtot;
Chris@42 63 static ptrdiff_t *local_ni=0, *local_starti=0;
Chris@42 64 static ptrdiff_t *local_no=0, *local_starto=0;
Chris@42 65 static ptrdiff_t *all_local_ni=0, *all_local_starti=0; /* n_pes x rnk arrays */
Chris@42 66 static ptrdiff_t *all_local_no=0, *all_local_starto=0; /* n_pes x rnk arrays */
Chris@42 67 static ptrdiff_t *istrides = 0, *ostrides = 0;
Chris@42 68 static ptrdiff_t *total_ni=0, *total_no=0;
Chris@42 69 static int *isend_cnt = 0, *isend_off = 0; /* for MPI_Scatterv */
Chris@42 70 static int *orecv_cnt = 0, *orecv_off = 0; /* for MPI_Gatherv */
Chris@42 71
Chris@42 72 static bench_real *local_in = 0, *local_out = 0;
Chris@42 73 static bench_real *all_local_in = 0, *all_local_out = 0;
Chris@42 74 static int all_local_in_alloc = 0, all_local_out_alloc = 0;
Chris@42 75 static FFTW(plan) plan_scramble_in = 0, plan_unscramble_out = 0;
Chris@42 76
Chris@42 77 static void alloc_rnk(int rnk_) {
Chris@42 78 rnk = rnk_;
Chris@42 79 bench_free(local_ni);
Chris@42 80 if (rnk == 0)
Chris@42 81 local_ni = 0;
Chris@42 82 else
Chris@42 83 local_ni = (ptrdiff_t *) bench_malloc(sizeof(ptrdiff_t) * rnk
Chris@42 84 * (8 + n_pes * 4));
Chris@42 85
Chris@42 86 local_starti = local_ni + rnk;
Chris@42 87 local_no = local_ni + 2 * rnk;
Chris@42 88 local_starto = local_ni + 3 * rnk;
Chris@42 89 istrides = local_ni + 4 * rnk;
Chris@42 90 ostrides = local_ni + 5 * rnk;
Chris@42 91 total_ni = local_ni + 6 * rnk;
Chris@42 92 total_no = local_ni + 7 * rnk;
Chris@42 93 all_local_ni = local_ni + 8 * rnk;
Chris@42 94 all_local_starti = local_ni + (8 + n_pes) * rnk;
Chris@42 95 all_local_no = local_ni + (8 + 2 * n_pes) * rnk;
Chris@42 96 all_local_starto = local_ni + (8 + 3 * n_pes) * rnk;
Chris@42 97 }
Chris@42 98
Chris@42 99 static void setup_gather_scatter(void)
Chris@42 100 {
Chris@42 101 int i, j;
Chris@42 102 ptrdiff_t off;
Chris@42 103
Chris@42 104 MPI_Gather(local_ni, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 105 all_local_ni, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 106 0, MPI_COMM_WORLD);
Chris@42 107 MPI_Bcast(all_local_ni, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
Chris@42 108 MPI_Gather(local_starti, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 109 all_local_starti, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 110 0, MPI_COMM_WORLD);
Chris@42 111 MPI_Bcast(all_local_starti, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
Chris@42 112
Chris@42 113 MPI_Gather(local_no, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 114 all_local_no, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 115 0, MPI_COMM_WORLD);
Chris@42 116 MPI_Bcast(all_local_no, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
Chris@42 117 MPI_Gather(local_starto, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 118 all_local_starto, rnk, FFTW_MPI_PTRDIFF_T,
Chris@42 119 0, MPI_COMM_WORLD);
Chris@42 120 MPI_Bcast(all_local_starto, rnk*n_pes, FFTW_MPI_PTRDIFF_T, 0, MPI_COMM_WORLD);
Chris@42 121
Chris@42 122 off = 0;
Chris@42 123 for (i = 0; i < n_pes; ++i) {
Chris@42 124 ptrdiff_t N = vn;
Chris@42 125 for (j = 0; j < rnk; ++j)
Chris@42 126 N *= all_local_ni[i * rnk + j];
Chris@42 127 isend_cnt[i] = N;
Chris@42 128 isend_off[i] = off;
Chris@42 129 off += N;
Chris@42 130 }
Chris@42 131 iNtot = off;
Chris@42 132 all_local_in_alloc = 1;
Chris@42 133
Chris@42 134 istrides[rnk - 1] = vn;
Chris@42 135 for (j = rnk - 2; j >= 0; --j)
Chris@42 136 istrides[j] = total_ni[j + 1] * istrides[j + 1];
Chris@42 137
Chris@42 138 off = 0;
Chris@42 139 for (i = 0; i < n_pes; ++i) {
Chris@42 140 ptrdiff_t N = vn;
Chris@42 141 for (j = 0; j < rnk; ++j)
Chris@42 142 N *= all_local_no[i * rnk + j];
Chris@42 143 orecv_cnt[i] = N;
Chris@42 144 orecv_off[i] = off;
Chris@42 145 off += N;
Chris@42 146 }
Chris@42 147 oNtot = off;
Chris@42 148 all_local_out_alloc = 1;
Chris@42 149
Chris@42 150 ostrides[rnk - 1] = vn;
Chris@42 151 for (j = rnk - 2; j >= 0; --j)
Chris@42 152 ostrides[j] = total_no[j + 1] * ostrides[j + 1];
Chris@42 153 }
Chris@42 154
Chris@42 155 static void copy_block_out(const bench_real *in,
Chris@42 156 int rnk, ptrdiff_t *n, ptrdiff_t *start,
Chris@42 157 ptrdiff_t is, ptrdiff_t *os, ptrdiff_t vn,
Chris@42 158 bench_real *out)
Chris@42 159 {
Chris@42 160 ptrdiff_t i;
Chris@42 161 if (rnk == 0) {
Chris@42 162 for (i = 0; i < vn; ++i)
Chris@42 163 out[i] = in[i];
Chris@42 164 }
Chris@42 165 else if (rnk == 1) { /* this case is just an optimization */
Chris@42 166 ptrdiff_t j;
Chris@42 167 out += start[0] * os[0];
Chris@42 168 for (j = 0; j < n[0]; ++j) {
Chris@42 169 for (i = 0; i < vn; ++i)
Chris@42 170 out[i] = in[i];
Chris@42 171 in += is;
Chris@42 172 out += os[0];
Chris@42 173 }
Chris@42 174 }
Chris@42 175 else {
Chris@42 176 /* we should do n[0] for locality, but this way is simpler to code */
Chris@42 177 for (i = 0; i < n[rnk - 1]; ++i)
Chris@42 178 copy_block_out(in + i * is,
Chris@42 179 rnk - 1, n, start, is * n[rnk - 1], os, vn,
Chris@42 180 out + (start[rnk - 1] + i) * os[rnk - 1]);
Chris@42 181 }
Chris@42 182 }
Chris@42 183
Chris@42 184 static void copy_block_in(bench_real *in,
Chris@42 185 int rnk, ptrdiff_t *n, ptrdiff_t *start,
Chris@42 186 ptrdiff_t is, ptrdiff_t *os, ptrdiff_t vn,
Chris@42 187 const bench_real *out)
Chris@42 188 {
Chris@42 189 ptrdiff_t i;
Chris@42 190 if (rnk == 0) {
Chris@42 191 for (i = 0; i < vn; ++i)
Chris@42 192 in[i] = out[i];
Chris@42 193 }
Chris@42 194 else if (rnk == 1) { /* this case is just an optimization */
Chris@42 195 ptrdiff_t j;
Chris@42 196 out += start[0] * os[0];
Chris@42 197 for (j = 0; j < n[0]; ++j) {
Chris@42 198 for (i = 0; i < vn; ++i)
Chris@42 199 in[i] = out[i];
Chris@42 200 in += is;
Chris@42 201 out += os[0];
Chris@42 202 }
Chris@42 203 }
Chris@42 204 else {
Chris@42 205 /* we should do n[0] for locality, but this way is simpler to code */
Chris@42 206 for (i = 0; i < n[rnk - 1]; ++i)
Chris@42 207 copy_block_in(in + i * is,
Chris@42 208 rnk - 1, n, start, is * n[rnk - 1], os, vn,
Chris@42 209 out + (start[rnk - 1] + i) * os[rnk - 1]);
Chris@42 210 }
Chris@42 211 }
Chris@42 212
Chris@42 213 static void do_scatter_in(bench_real *in)
Chris@42 214 {
Chris@42 215 bench_real *ali;
Chris@42 216 int i;
Chris@42 217 if (all_local_in_alloc) {
Chris@42 218 bench_free(all_local_in);
Chris@42 219 all_local_in = (bench_real*) bench_malloc(iNtot*sizeof(bench_real));
Chris@42 220 all_local_in_alloc = 0;
Chris@42 221 }
Chris@42 222 ali = all_local_in;
Chris@42 223 for (i = 0; i < n_pes; ++i) {
Chris@42 224 copy_block_in(ali,
Chris@42 225 rnk, all_local_ni + i * rnk,
Chris@42 226 all_local_starti + i * rnk,
Chris@42 227 vn, istrides, vn,
Chris@42 228 in);
Chris@42 229 ali += isend_cnt[i];
Chris@42 230 }
Chris@42 231 MPI_Scatterv(all_local_in, isend_cnt, isend_off, BENCH_MPI_TYPE,
Chris@42 232 local_in, isend_cnt[my_pe], BENCH_MPI_TYPE,
Chris@42 233 0, MPI_COMM_WORLD);
Chris@42 234 }
Chris@42 235
Chris@42 236 static void do_gather_out(bench_real *out)
Chris@42 237 {
Chris@42 238 bench_real *alo;
Chris@42 239 int i;
Chris@42 240
Chris@42 241 if (all_local_out_alloc) {
Chris@42 242 bench_free(all_local_out);
Chris@42 243 all_local_out = (bench_real*) bench_malloc(oNtot*sizeof(bench_real));
Chris@42 244 all_local_out_alloc = 0;
Chris@42 245 }
Chris@42 246 MPI_Gatherv(local_out, orecv_cnt[my_pe], BENCH_MPI_TYPE,
Chris@42 247 all_local_out, orecv_cnt, orecv_off, BENCH_MPI_TYPE,
Chris@42 248 0, MPI_COMM_WORLD);
Chris@42 249 MPI_Bcast(all_local_out, oNtot, BENCH_MPI_TYPE, 0, MPI_COMM_WORLD);
Chris@42 250 alo = all_local_out;
Chris@42 251 for (i = 0; i < n_pes; ++i) {
Chris@42 252 copy_block_out(alo,
Chris@42 253 rnk, all_local_no + i * rnk,
Chris@42 254 all_local_starto + i * rnk,
Chris@42 255 vn, ostrides, vn,
Chris@42 256 out);
Chris@42 257 alo += orecv_cnt[i];
Chris@42 258 }
Chris@42 259 }
Chris@42 260
Chris@42 261 static void alloc_local(ptrdiff_t nreal, int inplace)
Chris@42 262 {
Chris@42 263 bench_free(local_in);
Chris@42 264 if (local_out != local_in) bench_free(local_out);
Chris@42 265 local_in = local_out = 0;
Chris@42 266 if (nreal > 0) {
Chris@42 267 ptrdiff_t i;
Chris@42 268 local_in = (bench_real*) bench_malloc(nreal * sizeof(bench_real));
Chris@42 269 if (inplace)
Chris@42 270 local_out = local_in;
Chris@42 271 else
Chris@42 272 local_out = (bench_real*) bench_malloc(nreal * sizeof(bench_real));
Chris@42 273 for (i = 0; i < nreal; ++i) local_in[i] = local_out[i] = 0.0;
Chris@42 274 }
Chris@42 275 }
Chris@42 276
Chris@42 277 void after_problem_rcopy_from(bench_problem *p, bench_real *ri)
Chris@42 278 {
Chris@42 279 UNUSED(p);
Chris@42 280 do_scatter_in(ri);
Chris@42 281 if (plan_scramble_in) FFTW(execute)(plan_scramble_in);
Chris@42 282 }
Chris@42 283
Chris@42 284 void after_problem_rcopy_to(bench_problem *p, bench_real *ro)
Chris@42 285 {
Chris@42 286 UNUSED(p);
Chris@42 287 if (plan_unscramble_out) FFTW(execute)(plan_unscramble_out);
Chris@42 288 do_gather_out(ro);
Chris@42 289 }
Chris@42 290
Chris@42 291 void after_problem_ccopy_from(bench_problem *p, bench_real *ri, bench_real *ii)
Chris@42 292 {
Chris@42 293 UNUSED(ii);
Chris@42 294 after_problem_rcopy_from(p, ri);
Chris@42 295 }
Chris@42 296
Chris@42 297 void after_problem_ccopy_to(bench_problem *p, bench_real *ro, bench_real *io)
Chris@42 298 {
Chris@42 299 UNUSED(io);
Chris@42 300 after_problem_rcopy_to(p, ro);
Chris@42 301 }
Chris@42 302
Chris@42 303 void after_problem_hccopy_from(bench_problem *p, bench_real *ri, bench_real *ii)
Chris@42 304 {
Chris@42 305 UNUSED(ii);
Chris@42 306 after_problem_rcopy_from(p, ri);
Chris@42 307 }
Chris@42 308
Chris@42 309 void after_problem_hccopy_to(bench_problem *p, bench_real *ro, bench_real *io)
Chris@42 310 {
Chris@42 311 UNUSED(io);
Chris@42 312 after_problem_rcopy_to(p, ro);
Chris@42 313 }
Chris@42 314
Chris@42 315 static FFTW(plan) mkplan_transpose_local(ptrdiff_t nx, ptrdiff_t ny,
Chris@42 316 ptrdiff_t vn,
Chris@42 317 bench_real *in, bench_real *out)
Chris@42 318 {
Chris@42 319 FFTW(iodim64) hdims[3];
Chris@42 320 FFTW(r2r_kind) k[3];
Chris@42 321 FFTW(plan) pln;
Chris@42 322
Chris@42 323 hdims[0].n = nx;
Chris@42 324 hdims[0].is = ny * vn;
Chris@42 325 hdims[0].os = vn;
Chris@42 326 hdims[1].n = ny;
Chris@42 327 hdims[1].is = vn;
Chris@42 328 hdims[1].os = nx * vn;
Chris@42 329 hdims[2].n = vn;
Chris@42 330 hdims[2].is = 1;
Chris@42 331 hdims[2].os = 1;
Chris@42 332 k[0] = k[1] = k[2] = FFTW_R2HC;
Chris@42 333 pln = FFTW(plan_guru64_r2r)(0, 0, 3, hdims, in, out, k, FFTW_ESTIMATE);
Chris@42 334 BENCH_ASSERT(pln != 0);
Chris@42 335 return pln;
Chris@42 336 }
Chris@42 337
Chris@42 338 static int tensor_rowmajor_transposedp(bench_tensor *t)
Chris@42 339 {
Chris@42 340 bench_iodim *d;
Chris@42 341 int i;
Chris@42 342
Chris@42 343 BENCH_ASSERT(BENCH_FINITE_RNK(t->rnk));
Chris@42 344 if (t->rnk < 2)
Chris@42 345 return 0;
Chris@42 346
Chris@42 347 d = t->dims;
Chris@42 348 if (d[0].is != d[1].is * d[1].n
Chris@42 349 || d[0].os != d[1].is
Chris@42 350 || d[1].os != d[0].os * d[0].n)
Chris@42 351 return 0;
Chris@42 352 if (t->rnk > 2 && d[1].is != d[2].is * d[2].n)
Chris@42 353 return 0;
Chris@42 354 for (i = 2; i + 1 < t->rnk; ++i) {
Chris@42 355 d = t->dims + i;
Chris@42 356 if (d[0].is != d[1].is * d[1].n
Chris@42 357 || d[0].os != d[1].os * d[1].n)
Chris@42 358 return 0;
Chris@42 359 }
Chris@42 360
Chris@42 361 if (t->rnk > 2 && t->dims[t->rnk-1].is != t->dims[t->rnk-1].os)
Chris@42 362 return 0;
Chris@42 363 return 1;
Chris@42 364 }
Chris@42 365
Chris@42 366 static int tensor_contiguousp(bench_tensor *t, int s)
Chris@42 367 {
Chris@42 368 return (t->dims[t->rnk-1].is == s
Chris@42 369 && ((tensor_rowmajorp(t) &&
Chris@42 370 t->dims[t->rnk-1].is == t->dims[t->rnk-1].os)
Chris@42 371 || tensor_rowmajor_transposedp(t)));
Chris@42 372 }
Chris@42 373
Chris@42 374 static FFTW(plan) mkplan_complex(bench_problem *p, unsigned flags)
Chris@42 375 {
Chris@42 376 FFTW(plan) pln = 0;
Chris@42 377 int i;
Chris@42 378 ptrdiff_t ntot;
Chris@42 379
Chris@42 380 vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
Chris@42 381
Chris@42 382 if (p->sz->rnk < 1
Chris@42 383 || p->split
Chris@42 384 || !tensor_contiguousp(p->sz, vn)
Chris@42 385 || tensor_rowmajor_transposedp(p->sz)
Chris@42 386 || p->vecsz->rnk > 1
Chris@42 387 || (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
Chris@42 388 || p->vecsz->dims[0].os != 1)))
Chris@42 389 return 0;
Chris@42 390
Chris@42 391 alloc_rnk(p->sz->rnk);
Chris@42 392 for (i = 0; i < rnk; ++i) {
Chris@42 393 total_ni[i] = total_no[i] = p->sz->dims[i].n;
Chris@42 394 local_ni[i] = local_no[i] = total_ni[i];
Chris@42 395 local_starti[i] = local_starto[i] = 0;
Chris@42 396 }
Chris@42 397 if (rnk > 1) {
Chris@42 398 ptrdiff_t n, start, nT, startT;
Chris@42 399 ntot = FFTW(mpi_local_size_many_transposed)
Chris@42 400 (p->sz->rnk, total_ni, vn,
Chris@42 401 FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
Chris@42 402 MPI_COMM_WORLD,
Chris@42 403 &n, &start, &nT, &startT);
Chris@42 404 if (flags & FFTW_MPI_TRANSPOSED_IN) {
Chris@42 405 local_ni[1] = nT;
Chris@42 406 local_starti[1] = startT;
Chris@42 407 }
Chris@42 408 else {
Chris@42 409 local_ni[0] = n;
Chris@42 410 local_starti[0] = start;
Chris@42 411 }
Chris@42 412 if (flags & FFTW_MPI_TRANSPOSED_OUT) {
Chris@42 413 local_no[1] = nT;
Chris@42 414 local_starto[1] = startT;
Chris@42 415 }
Chris@42 416 else {
Chris@42 417 local_no[0] = n;
Chris@42 418 local_starto[0] = start;
Chris@42 419 }
Chris@42 420 }
Chris@42 421 else if (rnk == 1) {
Chris@42 422 ntot = FFTW(mpi_local_size_many_1d)
Chris@42 423 (total_ni[0], vn, MPI_COMM_WORLD, p->sign, flags,
Chris@42 424 local_ni, local_starti, local_no, local_starto);
Chris@42 425 }
Chris@42 426 alloc_local(ntot * 2, p->in == p->out);
Chris@42 427
Chris@42 428 pln = FFTW(mpi_plan_many_dft)(p->sz->rnk, total_ni, vn,
Chris@42 429 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 430 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 431 (FFTW(complex) *) local_in,
Chris@42 432 (FFTW(complex) *) local_out,
Chris@42 433 MPI_COMM_WORLD, p->sign, flags);
Chris@42 434
Chris@42 435 vn *= 2;
Chris@42 436
Chris@42 437 if (rnk > 1) {
Chris@42 438 ptrdiff_t nrest = 1;
Chris@42 439 for (i = 2; i < rnk; ++i) nrest *= p->sz->dims[i].n;
Chris@42 440 if (flags & FFTW_MPI_TRANSPOSED_IN)
Chris@42 441 plan_scramble_in = mkplan_transpose_local(
Chris@42 442 p->sz->dims[0].n, local_ni[1], vn * nrest,
Chris@42 443 local_in, local_in);
Chris@42 444 if (flags & FFTW_MPI_TRANSPOSED_OUT)
Chris@42 445 plan_unscramble_out = mkplan_transpose_local(
Chris@42 446 local_no[1], p->sz->dims[0].n, vn * nrest,
Chris@42 447 local_out, local_out);
Chris@42 448 }
Chris@42 449
Chris@42 450 return pln;
Chris@42 451 }
Chris@42 452
Chris@42 453 static int tensor_real_contiguousp(bench_tensor *t, int sign, int s)
Chris@42 454 {
Chris@42 455 return (t->dims[t->rnk-1].is == s
Chris@42 456 && ((tensor_real_rowmajorp(t, sign, 1) &&
Chris@42 457 t->dims[t->rnk-1].is == t->dims[t->rnk-1].os)));
Chris@42 458 }
Chris@42 459
Chris@42 460 static FFTW(plan) mkplan_real(bench_problem *p, unsigned flags)
Chris@42 461 {
Chris@42 462 FFTW(plan) pln = 0;
Chris@42 463 int i;
Chris@42 464 ptrdiff_t ntot;
Chris@42 465
Chris@42 466 vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
Chris@42 467
Chris@42 468 if (p->sz->rnk < 2
Chris@42 469 || p->split
Chris@42 470 || !tensor_real_contiguousp(p->sz, p->sign, vn)
Chris@42 471 || tensor_rowmajor_transposedp(p->sz)
Chris@42 472 || p->vecsz->rnk > 1
Chris@42 473 || (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
Chris@42 474 || p->vecsz->dims[0].os != 1)))
Chris@42 475 return 0;
Chris@42 476
Chris@42 477 alloc_rnk(p->sz->rnk);
Chris@42 478 for (i = 0; i < rnk; ++i) {
Chris@42 479 total_ni[i] = total_no[i] = p->sz->dims[i].n;
Chris@42 480 local_ni[i] = local_no[i] = total_ni[i];
Chris@42 481 local_starti[i] = local_starto[i] = 0;
Chris@42 482 }
Chris@42 483 local_ni[rnk-1] = local_no[rnk-1] = total_ni[rnk-1] = total_no[rnk-1]
Chris@42 484 = p->sz->dims[rnk-1].n / 2 + 1;
Chris@42 485 {
Chris@42 486 ptrdiff_t n, start, nT, startT;
Chris@42 487 ntot = FFTW(mpi_local_size_many_transposed)
Chris@42 488 (p->sz->rnk, total_ni, vn,
Chris@42 489 FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
Chris@42 490 MPI_COMM_WORLD,
Chris@42 491 &n, &start, &nT, &startT);
Chris@42 492 if (flags & FFTW_MPI_TRANSPOSED_IN) {
Chris@42 493 local_ni[1] = nT;
Chris@42 494 local_starti[1] = startT;
Chris@42 495 }
Chris@42 496 else {
Chris@42 497 local_ni[0] = n;
Chris@42 498 local_starti[0] = start;
Chris@42 499 }
Chris@42 500 if (flags & FFTW_MPI_TRANSPOSED_OUT) {
Chris@42 501 local_no[1] = nT;
Chris@42 502 local_starto[1] = startT;
Chris@42 503 }
Chris@42 504 else {
Chris@42 505 local_no[0] = n;
Chris@42 506 local_starto[0] = start;
Chris@42 507 }
Chris@42 508 }
Chris@42 509 alloc_local(ntot * 2, p->in == p->out);
Chris@42 510
Chris@42 511 total_ni[rnk - 1] = p->sz->dims[rnk - 1].n;
Chris@42 512 if (p->sign < 0)
Chris@42 513 pln = FFTW(mpi_plan_many_dft_r2c)(p->sz->rnk, total_ni, vn,
Chris@42 514 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 515 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 516 local_in,
Chris@42 517 (FFTW(complex) *) local_out,
Chris@42 518 MPI_COMM_WORLD, flags);
Chris@42 519 else
Chris@42 520 pln = FFTW(mpi_plan_many_dft_c2r)(p->sz->rnk, total_ni, vn,
Chris@42 521 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 522 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 523 (FFTW(complex) *) local_in,
Chris@42 524 local_out,
Chris@42 525 MPI_COMM_WORLD, flags);
Chris@42 526
Chris@42 527 total_ni[rnk - 1] = p->sz->dims[rnk - 1].n / 2 + 1;
Chris@42 528 vn *= 2;
Chris@42 529
Chris@42 530 {
Chris@42 531 ptrdiff_t nrest = 1;
Chris@42 532 for (i = 2; i < rnk; ++i) nrest *= total_ni[i];
Chris@42 533 if (flags & FFTW_MPI_TRANSPOSED_IN)
Chris@42 534 plan_scramble_in = mkplan_transpose_local(
Chris@42 535 total_ni[0], local_ni[1], vn * nrest,
Chris@42 536 local_in, local_in);
Chris@42 537 if (flags & FFTW_MPI_TRANSPOSED_OUT)
Chris@42 538 plan_unscramble_out = mkplan_transpose_local(
Chris@42 539 local_no[1], total_ni[0], vn * nrest,
Chris@42 540 local_out, local_out);
Chris@42 541 }
Chris@42 542
Chris@42 543 return pln;
Chris@42 544 }
Chris@42 545
Chris@42 546 static FFTW(plan) mkplan_transpose(bench_problem *p, unsigned flags)
Chris@42 547 {
Chris@42 548 ptrdiff_t ntot, nx, ny;
Chris@42 549 int ix=0, iy=1, i;
Chris@42 550 const bench_iodim *d = p->vecsz->dims;
Chris@42 551 FFTW(plan) pln;
Chris@42 552
Chris@42 553 if (p->vecsz->rnk == 3) {
Chris@42 554 for (i = 0; i < 3; ++i)
Chris@42 555 if (d[i].is == 1 && d[i].os == 1) {
Chris@42 556 vn = d[i].n;
Chris@42 557 ix = (i + 1) % 3;
Chris@42 558 iy = (i + 2) % 3;
Chris@42 559 break;
Chris@42 560 }
Chris@42 561 if (i == 3) return 0;
Chris@42 562 }
Chris@42 563 else {
Chris@42 564 vn = 1;
Chris@42 565 ix = 0;
Chris@42 566 iy = 1;
Chris@42 567 }
Chris@42 568
Chris@42 569 if (d[ix].is == d[iy].n * vn && d[ix].os == vn
Chris@42 570 && d[iy].os == d[ix].n * vn && d[iy].is == vn) {
Chris@42 571 nx = d[ix].n;
Chris@42 572 ny = d[iy].n;
Chris@42 573 }
Chris@42 574 else if (d[iy].is == d[ix].n * vn && d[iy].os == vn
Chris@42 575 && d[ix].os == d[iy].n * vn && d[ix].is == vn) {
Chris@42 576 nx = d[iy].n;
Chris@42 577 ny = d[ix].n;
Chris@42 578 }
Chris@42 579 else
Chris@42 580 return 0;
Chris@42 581
Chris@42 582 alloc_rnk(2);
Chris@42 583 ntot = vn * FFTW(mpi_local_size_2d_transposed)(nx, ny, MPI_COMM_WORLD,
Chris@42 584 &local_ni[0],
Chris@42 585 &local_starti[0],
Chris@42 586 &local_no[0],
Chris@42 587 &local_starto[0]);
Chris@42 588 local_ni[1] = ny;
Chris@42 589 local_starti[1] = 0;
Chris@42 590 local_no[1] = nx;
Chris@42 591 local_starto[1] = 0;
Chris@42 592 total_ni[0] = nx; total_ni[1] = ny;
Chris@42 593 total_no[1] = nx; total_no[0] = ny;
Chris@42 594 alloc_local(ntot, p->in == p->out);
Chris@42 595
Chris@42 596 pln = FFTW(mpi_plan_many_transpose)(nx, ny, vn,
Chris@42 597 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 598 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 599 local_in, local_out,
Chris@42 600 MPI_COMM_WORLD, flags);
Chris@42 601
Chris@42 602 if (flags & FFTW_MPI_TRANSPOSED_IN)
Chris@42 603 plan_scramble_in = mkplan_transpose_local(local_ni[0], ny, vn,
Chris@42 604 local_in, local_in);
Chris@42 605 if (flags & FFTW_MPI_TRANSPOSED_OUT)
Chris@42 606 plan_unscramble_out = mkplan_transpose_local
Chris@42 607 (nx, local_no[0], vn, local_out, local_out);
Chris@42 608
Chris@42 609 #if 0
Chris@42 610 if (pln && vn == 1) {
Chris@42 611 int i, j;
Chris@42 612 bench_real *ri = (bench_real *) p->in;
Chris@42 613 bench_real *ro = (bench_real *) p->out;
Chris@42 614 if (!ri || !ro) return pln;
Chris@42 615 setup_gather_scatter();
Chris@42 616 for (i = 0; i < nx * ny; ++i)
Chris@42 617 ri[i] = i;
Chris@42 618 after_problem_rcopy_from(p, ri);
Chris@42 619 FFTW(execute)(pln);
Chris@42 620 after_problem_rcopy_to(p, ro);
Chris@42 621 if (my_pe == 0) {
Chris@42 622 for (i = 0; i < nx; ++i) {
Chris@42 623 for (j = 0; j < ny; ++j)
Chris@42 624 printf(" %3g", ro[j * nx + i]);
Chris@42 625 printf("\n");
Chris@42 626 }
Chris@42 627 }
Chris@42 628 }
Chris@42 629 #endif
Chris@42 630
Chris@42 631 return pln;
Chris@42 632 }
Chris@42 633
Chris@42 634 static FFTW(plan) mkplan_r2r(bench_problem *p, unsigned flags)
Chris@42 635 {
Chris@42 636 FFTW(plan) pln = 0;
Chris@42 637 int i;
Chris@42 638 ptrdiff_t ntot;
Chris@42 639 FFTW(r2r_kind) *k;
Chris@42 640
Chris@42 641 if ((p->sz->rnk == 0 || (p->sz->rnk == 1 && p->sz->dims[0].n == 1))
Chris@42 642 && p->vecsz->rnk >= 2 && p->vecsz->rnk <= 3)
Chris@42 643 return mkplan_transpose(p, flags);
Chris@42 644
Chris@42 645 vn = p->vecsz->rnk == 1 ? p->vecsz->dims[0].n : 1;
Chris@42 646
Chris@42 647 if (p->sz->rnk < 1
Chris@42 648 || p->split
Chris@42 649 || !tensor_contiguousp(p->sz, vn)
Chris@42 650 || tensor_rowmajor_transposedp(p->sz)
Chris@42 651 || p->vecsz->rnk > 1
Chris@42 652 || (p->vecsz->rnk == 1 && (p->vecsz->dims[0].is != 1
Chris@42 653 || p->vecsz->dims[0].os != 1)))
Chris@42 654 return 0;
Chris@42 655
Chris@42 656 alloc_rnk(p->sz->rnk);
Chris@42 657 for (i = 0; i < rnk; ++i) {
Chris@42 658 total_ni[i] = total_no[i] = p->sz->dims[i].n;
Chris@42 659 local_ni[i] = local_no[i] = total_ni[i];
Chris@42 660 local_starti[i] = local_starto[i] = 0;
Chris@42 661 }
Chris@42 662 if (rnk > 1) {
Chris@42 663 ptrdiff_t n, start, nT, startT;
Chris@42 664 ntot = FFTW(mpi_local_size_many_transposed)
Chris@42 665 (p->sz->rnk, total_ni, vn,
Chris@42 666 FFTW_MPI_DEFAULT_BLOCK, FFTW_MPI_DEFAULT_BLOCK,
Chris@42 667 MPI_COMM_WORLD,
Chris@42 668 &n, &start, &nT, &startT);
Chris@42 669 if (flags & FFTW_MPI_TRANSPOSED_IN) {
Chris@42 670 local_ni[1] = nT;
Chris@42 671 local_starti[1] = startT;
Chris@42 672 }
Chris@42 673 else {
Chris@42 674 local_ni[0] = n;
Chris@42 675 local_starti[0] = start;
Chris@42 676 }
Chris@42 677 if (flags & FFTW_MPI_TRANSPOSED_OUT) {
Chris@42 678 local_no[1] = nT;
Chris@42 679 local_starto[1] = startT;
Chris@42 680 }
Chris@42 681 else {
Chris@42 682 local_no[0] = n;
Chris@42 683 local_starto[0] = start;
Chris@42 684 }
Chris@42 685 }
Chris@42 686 else if (rnk == 1) {
Chris@42 687 ntot = FFTW(mpi_local_size_many_1d)
Chris@42 688 (total_ni[0], vn, MPI_COMM_WORLD, p->sign, flags,
Chris@42 689 local_ni, local_starti, local_no, local_starto);
Chris@42 690 }
Chris@42 691 alloc_local(ntot, p->in == p->out);
Chris@42 692
Chris@42 693 k = (FFTW(r2r_kind) *) bench_malloc(sizeof(FFTW(r2r_kind)) * p->sz->rnk);
Chris@42 694 for (i = 0; i < p->sz->rnk; ++i)
Chris@42 695 switch (p->k[i]) {
Chris@42 696 case R2R_R2HC: k[i] = FFTW_R2HC; break;
Chris@42 697 case R2R_HC2R: k[i] = FFTW_HC2R; break;
Chris@42 698 case R2R_DHT: k[i] = FFTW_DHT; break;
Chris@42 699 case R2R_REDFT00: k[i] = FFTW_REDFT00; break;
Chris@42 700 case R2R_REDFT01: k[i] = FFTW_REDFT01; break;
Chris@42 701 case R2R_REDFT10: k[i] = FFTW_REDFT10; break;
Chris@42 702 case R2R_REDFT11: k[i] = FFTW_REDFT11; break;
Chris@42 703 case R2R_RODFT00: k[i] = FFTW_RODFT00; break;
Chris@42 704 case R2R_RODFT01: k[i] = FFTW_RODFT01; break;
Chris@42 705 case R2R_RODFT10: k[i] = FFTW_RODFT10; break;
Chris@42 706 case R2R_RODFT11: k[i] = FFTW_RODFT11; break;
Chris@42 707 default: BENCH_ASSERT(0);
Chris@42 708 }
Chris@42 709
Chris@42 710 pln = FFTW(mpi_plan_many_r2r)(p->sz->rnk, total_ni, vn,
Chris@42 711 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 712 FFTW_MPI_DEFAULT_BLOCK,
Chris@42 713 local_in, local_out,
Chris@42 714 MPI_COMM_WORLD, k, flags);
Chris@42 715 bench_free(k);
Chris@42 716
Chris@42 717 if (rnk > 1) {
Chris@42 718 ptrdiff_t nrest = 1;
Chris@42 719 for (i = 2; i < rnk; ++i) nrest *= p->sz->dims[i].n;
Chris@42 720 if (flags & FFTW_MPI_TRANSPOSED_IN)
Chris@42 721 plan_scramble_in = mkplan_transpose_local(
Chris@42 722 p->sz->dims[0].n, local_ni[1], vn * nrest,
Chris@42 723 local_in, local_in);
Chris@42 724 if (flags & FFTW_MPI_TRANSPOSED_OUT)
Chris@42 725 plan_unscramble_out = mkplan_transpose_local(
Chris@42 726 local_no[1], p->sz->dims[0].n, vn * nrest,
Chris@42 727 local_out, local_out);
Chris@42 728 }
Chris@42 729
Chris@42 730 return pln;
Chris@42 731 }
Chris@42 732
Chris@42 733 FFTW(plan) mkplan(bench_problem *p, unsigned flags)
Chris@42 734 {
Chris@42 735 FFTW(plan) pln = 0;
Chris@42 736 FFTW(destroy_plan)(plan_scramble_in); plan_scramble_in = 0;
Chris@42 737 FFTW(destroy_plan)(plan_unscramble_out); plan_unscramble_out = 0;
Chris@42 738 if (p->scrambled_in) {
Chris@42 739 if (p->sz->rnk == 1 && p->sz->dims[0].n != 1)
Chris@42 740 flags |= FFTW_MPI_SCRAMBLED_IN;
Chris@42 741 else
Chris@42 742 flags |= FFTW_MPI_TRANSPOSED_IN;
Chris@42 743 }
Chris@42 744 if (p->scrambled_out) {
Chris@42 745 if (p->sz->rnk == 1 && p->sz->dims[0].n != 1)
Chris@42 746 flags |= FFTW_MPI_SCRAMBLED_OUT;
Chris@42 747 else
Chris@42 748 flags |= FFTW_MPI_TRANSPOSED_OUT;
Chris@42 749 }
Chris@42 750 switch (p->kind) {
Chris@42 751 case PROBLEM_COMPLEX:
Chris@42 752 pln =mkplan_complex(p, flags);
Chris@42 753 break;
Chris@42 754 case PROBLEM_REAL:
Chris@42 755 pln = mkplan_real(p, flags);
Chris@42 756 break;
Chris@42 757 case PROBLEM_R2R:
Chris@42 758 pln = mkplan_r2r(p, flags);
Chris@42 759 break;
Chris@42 760 default: BENCH_ASSERT(0);
Chris@42 761 }
Chris@42 762 if (pln) setup_gather_scatter();
Chris@42 763 return pln;
Chris@42 764 }
Chris@42 765
Chris@42 766 void main_init(int *argc, char ***argv)
Chris@42 767 {
Chris@42 768 #ifdef HAVE_SMP
Chris@42 769 # if MPI_VERSION >= 2 /* for MPI_Init_thread */
Chris@42 770 int provided;
Chris@42 771 MPI_Init_thread(argc, argv, MPI_THREAD_FUNNELED, &provided);
Chris@42 772 threads_ok = provided >= MPI_THREAD_FUNNELED;
Chris@42 773 # else
Chris@42 774 MPI_Init(argc, argv);
Chris@42 775 threads_ok = 0;
Chris@42 776 # endif
Chris@42 777 #else
Chris@42 778 MPI_Init(argc, argv);
Chris@42 779 #endif
Chris@42 780 MPI_Comm_rank(MPI_COMM_WORLD, &my_pe);
Chris@42 781 MPI_Comm_size(MPI_COMM_WORLD, &n_pes);
Chris@42 782 if (my_pe != 0) verbose = -999;
Chris@42 783 no_speed_allocation = 1; /* so we can benchmark transforms > memory */
Chris@42 784 always_pad_real = 1; /* out-of-place real transforms are padded */
Chris@42 785 isend_cnt = (int *) bench_malloc(sizeof(int) * n_pes);
Chris@42 786 isend_off = (int *) bench_malloc(sizeof(int) * n_pes);
Chris@42 787 orecv_cnt = (int *) bench_malloc(sizeof(int) * n_pes);
Chris@42 788 orecv_off = (int *) bench_malloc(sizeof(int) * n_pes);
Chris@42 789
Chris@42 790 /* init_threads must be called before any other FFTW function,
Chris@42 791 including mpi_init, because it has to register the threads hooks
Chris@42 792 before the planner is initalized */
Chris@42 793 #ifdef HAVE_SMP
Chris@42 794 if (threads_ok) { BENCH_ASSERT(FFTW(init_threads)()); }
Chris@42 795 #endif
Chris@42 796 FFTW(mpi_init)();
Chris@42 797 }
Chris@42 798
Chris@42 799 void initial_cleanup(void)
Chris@42 800 {
Chris@42 801 alloc_rnk(0);
Chris@42 802 alloc_local(0, 0);
Chris@42 803 bench_free(all_local_in); all_local_in = 0;
Chris@42 804 bench_free(all_local_out); all_local_out = 0;
Chris@42 805 bench_free(isend_off); isend_off = 0;
Chris@42 806 bench_free(isend_cnt); isend_cnt = 0;
Chris@42 807 bench_free(orecv_off); orecv_off = 0;
Chris@42 808 bench_free(orecv_cnt); orecv_cnt = 0;
Chris@42 809 FFTW(destroy_plan)(plan_scramble_in); plan_scramble_in = 0;
Chris@42 810 FFTW(destroy_plan)(plan_unscramble_out); plan_unscramble_out = 0;
Chris@42 811 }
Chris@42 812
Chris@42 813 void final_cleanup(void)
Chris@42 814 {
Chris@42 815 MPI_Finalize();
Chris@42 816 }
Chris@42 817
Chris@42 818 void bench_exit(int status)
Chris@42 819 {
Chris@42 820 MPI_Abort(MPI_COMM_WORLD, status);
Chris@42 821 }
Chris@42 822
Chris@42 823 double bench_cost_postprocess(double cost)
Chris@42 824 {
Chris@42 825 double cost_max;
Chris@42 826 MPI_Allreduce(&cost, &cost_max, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
Chris@42 827 return cost_max;
Chris@42 828 }
Chris@42 829
Chris@42 830
Chris@42 831 int import_wisdom(FILE *f)
Chris@42 832 {
Chris@42 833 int success = 1, sall;
Chris@42 834 if (my_pe == 0) success = FFTW(import_wisdom_from_file)(f);
Chris@42 835 FFTW(mpi_broadcast_wisdom)(MPI_COMM_WORLD);
Chris@42 836 MPI_Allreduce(&success, &sall, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
Chris@42 837 return sall;
Chris@42 838 }
Chris@42 839
Chris@42 840 void export_wisdom(FILE *f)
Chris@42 841 {
Chris@42 842 FFTW(mpi_gather_wisdom)(MPI_COMM_WORLD);
Chris@42 843 if (my_pe == 0) FFTW(export_wisdom_to_file)(f);
Chris@42 844 }