annotate fft/fftw/fftw-3.3.4/mpi/mpi-bench.c @ 40:223f770b5341 kissfft-double tip

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