annotate src/fftw-3.3.3/mpi/mpi-bench.c @ 169:223a55898ab9 tip default

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