annotate src/fftw-3.3.3/mpi/api.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 * Copyright (c) 2003, 2007-11 Matteo Frigo
cannam@95 3 * Copyright (c) 2003, 2007-11 Massachusetts Institute of Technology
cannam@95 4 *
cannam@95 5 * This program is free software; you can redistribute it and/or modify
cannam@95 6 * it under the terms of the GNU General Public License as published by
cannam@95 7 * the Free Software Foundation; either version 2 of the License, or
cannam@95 8 * (at your option) any later version.
cannam@95 9 *
cannam@95 10 * This program is distributed in the hope that it will be useful,
cannam@95 11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
cannam@95 12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
cannam@95 13 * GNU General Public License for more details.
cannam@95 14 *
cannam@95 15 * You should have received a copy of the GNU General Public License
cannam@95 16 * along with this program; if not, write to the Free Software
cannam@95 17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
cannam@95 18 *
cannam@95 19 */
cannam@95 20
cannam@95 21 #include "api.h"
cannam@95 22 #include "fftw3-mpi.h"
cannam@95 23 #include "ifftw-mpi.h"
cannam@95 24 #include "mpi-transpose.h"
cannam@95 25 #include "mpi-dft.h"
cannam@95 26 #include "mpi-rdft.h"
cannam@95 27 #include "mpi-rdft2.h"
cannam@95 28
cannam@95 29 /* Convert API flags to internal MPI flags. */
cannam@95 30 #define MPI_FLAGS(f) ((f) >> 27)
cannam@95 31
cannam@95 32 /*************************************************************************/
cannam@95 33
cannam@95 34 static int mpi_inited = 0;
cannam@95 35
cannam@95 36 static MPI_Comm problem_comm(const problem *p) {
cannam@95 37 switch (p->adt->problem_kind) {
cannam@95 38 case PROBLEM_MPI_DFT:
cannam@95 39 return ((const problem_mpi_dft *) p)->comm;
cannam@95 40 case PROBLEM_MPI_RDFT:
cannam@95 41 return ((const problem_mpi_rdft *) p)->comm;
cannam@95 42 case PROBLEM_MPI_RDFT2:
cannam@95 43 return ((const problem_mpi_rdft2 *) p)->comm;
cannam@95 44 case PROBLEM_MPI_TRANSPOSE:
cannam@95 45 return ((const problem_mpi_transpose *) p)->comm;
cannam@95 46 default:
cannam@95 47 return MPI_COMM_NULL;
cannam@95 48 }
cannam@95 49 }
cannam@95 50
cannam@95 51 /* used to synchronize cost measurements (timing or estimation)
cannam@95 52 across all processes for an MPI problem, which is critical to
cannam@95 53 ensure that all processes decide to use the same MPI plans
cannam@95 54 (whereas serial plans need not be syncronized). */
cannam@95 55 static double cost_hook(const problem *p, double t, cost_kind k)
cannam@95 56 {
cannam@95 57 MPI_Comm comm = problem_comm(p);
cannam@95 58 double tsum;
cannam@95 59 if (comm == MPI_COMM_NULL) return t;
cannam@95 60 MPI_Allreduce(&t, &tsum, 1, MPI_DOUBLE,
cannam@95 61 k == COST_SUM ? MPI_SUM : MPI_MAX, comm);
cannam@95 62 return tsum;
cannam@95 63 }
cannam@95 64
cannam@95 65 /* Used to reject wisdom that is not in sync across all processes
cannam@95 66 for an MPI problem, which is critical to ensure that all processes
cannam@95 67 decide to use the same MPI plans. (Even though costs are synchronized,
cannam@95 68 above, out-of-sync wisdom may result from plans being produced
cannam@95 69 by communicators that do not span all processes, either from a
cannam@95 70 user-specified communicator or e.g. from transpose-recurse. */
cannam@95 71 static int wisdom_ok_hook(const problem *p, flags_t flags)
cannam@95 72 {
cannam@95 73 MPI_Comm comm = problem_comm(p);
cannam@95 74 int eq_me, eq_all;
cannam@95 75 /* unpack flags bitfield, since MPI communications may involve
cannam@95 76 byte-order changes and MPI cannot do this for bit fields */
cannam@95 77 #if SIZEOF_UNSIGNED_INT >= 4 /* must be big enough to hold 20-bit fields */
cannam@95 78 unsigned int f[5];
cannam@95 79 #else
cannam@95 80 unsigned long f[5]; /* at least 32 bits as per C standard */
cannam@95 81 #endif
cannam@95 82
cannam@95 83 if (comm == MPI_COMM_NULL) return 1; /* non-MPI wisdom is always ok */
cannam@95 84
cannam@95 85 if (XM(any_true)(0, comm)) return 0; /* some process had nowisdom_hook */
cannam@95 86
cannam@95 87 /* otherwise, check that the flags and solver index are identical
cannam@95 88 on all processes in this problem's communicator.
cannam@95 89
cannam@95 90 TO DO: possibly we can relax strict equality, but it is
cannam@95 91 critical to ensure that any flags which affect what plan is
cannam@95 92 created (and whether the solver is applicable) are the same,
cannam@95 93 e.g. DESTROY_INPUT, NO_UGLY, etcetera. (If the MPI algorithm
cannam@95 94 differs between processes, deadlocks/crashes generally result.) */
cannam@95 95 f[0] = flags.l;
cannam@95 96 f[1] = flags.hash_info;
cannam@95 97 f[2] = flags.timelimit_impatience;
cannam@95 98 f[3] = flags.u;
cannam@95 99 f[4] = flags.slvndx;
cannam@95 100 MPI_Bcast(f, 5,
cannam@95 101 SIZEOF_UNSIGNED_INT >= 4 ? MPI_UNSIGNED : MPI_UNSIGNED_LONG,
cannam@95 102 0, comm);
cannam@95 103 eq_me = f[0] == flags.l && f[1] == flags.hash_info
cannam@95 104 && f[2] == flags.timelimit_impatience
cannam@95 105 && f[3] == flags.u && f[4] == flags.slvndx;
cannam@95 106 MPI_Allreduce(&eq_me, &eq_all, 1, MPI_INT, MPI_LAND, comm);
cannam@95 107 return eq_all;
cannam@95 108 }
cannam@95 109
cannam@95 110 /* This hook is called when wisdom is not found. The any_true here
cannam@95 111 matches up with the any_true in wisdom_ok_hook, in order to handle
cannam@95 112 the case where some processes had wisdom (and called wisdom_ok_hook)
cannam@95 113 and some processes didn't have wisdom (and called nowisdom_hook). */
cannam@95 114 static void nowisdom_hook(const problem *p)
cannam@95 115 {
cannam@95 116 MPI_Comm comm = problem_comm(p);
cannam@95 117 if (comm == MPI_COMM_NULL) return; /* nothing to do for non-MPI p */
cannam@95 118 XM(any_true)(1, comm); /* signal nowisdom to any wisdom_ok_hook */
cannam@95 119 }
cannam@95 120
cannam@95 121 /* needed to synchronize planner bogosity flag, in case non-MPI problems
cannam@95 122 on a subset of processes encountered bogus wisdom */
cannam@95 123 static wisdom_state_t bogosity_hook(wisdom_state_t state, const problem *p)
cannam@95 124 {
cannam@95 125 MPI_Comm comm = problem_comm(p);
cannam@95 126 if (comm != MPI_COMM_NULL /* an MPI problem */
cannam@95 127 && XM(any_true)(state == WISDOM_IS_BOGUS, comm)) /* bogus somewhere */
cannam@95 128 return WISDOM_IS_BOGUS;
cannam@95 129 return state;
cannam@95 130 }
cannam@95 131
cannam@95 132 void XM(init)(void)
cannam@95 133 {
cannam@95 134 if (!mpi_inited) {
cannam@95 135 planner *plnr = X(the_planner)();
cannam@95 136 plnr->cost_hook = cost_hook;
cannam@95 137 plnr->wisdom_ok_hook = wisdom_ok_hook;
cannam@95 138 plnr->nowisdom_hook = nowisdom_hook;
cannam@95 139 plnr->bogosity_hook = bogosity_hook;
cannam@95 140 XM(conf_standard)(plnr);
cannam@95 141 mpi_inited = 1;
cannam@95 142 }
cannam@95 143 }
cannam@95 144
cannam@95 145 void XM(cleanup)(void)
cannam@95 146 {
cannam@95 147 X(cleanup)();
cannam@95 148 mpi_inited = 0;
cannam@95 149 }
cannam@95 150
cannam@95 151 /*************************************************************************/
cannam@95 152
cannam@95 153 static dtensor *mkdtensor_api(int rnk, const XM(ddim) *dims0)
cannam@95 154 {
cannam@95 155 dtensor *x = XM(mkdtensor)(rnk);
cannam@95 156 int i;
cannam@95 157 for (i = 0; i < rnk; ++i) {
cannam@95 158 x->dims[i].n = dims0[i].n;
cannam@95 159 x->dims[i].b[IB] = dims0[i].ib;
cannam@95 160 x->dims[i].b[OB] = dims0[i].ob;
cannam@95 161 }
cannam@95 162 return x;
cannam@95 163 }
cannam@95 164
cannam@95 165 static dtensor *default_sz(int rnk, const XM(ddim) *dims0, int n_pes,
cannam@95 166 int rdft2)
cannam@95 167 {
cannam@95 168 dtensor *sz = XM(mkdtensor)(rnk);
cannam@95 169 dtensor *sz0 = mkdtensor_api(rnk, dims0);
cannam@95 170 block_kind k;
cannam@95 171 int i;
cannam@95 172
cannam@95 173 for (i = 0; i < rnk; ++i)
cannam@95 174 sz->dims[i].n = dims0[i].n;
cannam@95 175
cannam@95 176 if (rdft2) sz->dims[rnk-1].n = dims0[rnk-1].n / 2 + 1;
cannam@95 177
cannam@95 178 for (i = 0; i < rnk; ++i) {
cannam@95 179 sz->dims[i].b[IB] = dims0[i].ib ? dims0[i].ib : sz->dims[i].n;
cannam@95 180 sz->dims[i].b[OB] = dims0[i].ob ? dims0[i].ob : sz->dims[i].n;
cannam@95 181 }
cannam@95 182
cannam@95 183 /* If we haven't used all of the processes yet, and some of the
cannam@95 184 block sizes weren't specified (i.e. 0), then set the
cannam@95 185 unspecified blocks so as to use as many processes as
cannam@95 186 possible with as few distributed dimensions as possible. */
cannam@95 187 FORALL_BLOCK_KIND(k) {
cannam@95 188 INT nb = XM(num_blocks_total)(sz, k);
cannam@95 189 INT np = n_pes / nb;
cannam@95 190 for (i = 0; i < rnk && np > 1; ++i)
cannam@95 191 if (!sz0->dims[i].b[k]) {
cannam@95 192 sz->dims[i].b[k] = XM(default_block)(sz->dims[i].n, np);
cannam@95 193 nb *= XM(num_blocks)(sz->dims[i].n, sz->dims[i].b[k]);
cannam@95 194 np = n_pes / nb;
cannam@95 195 }
cannam@95 196 }
cannam@95 197
cannam@95 198 if (rdft2) sz->dims[rnk-1].n = dims0[rnk-1].n;
cannam@95 199
cannam@95 200 /* punt for 1d prime */
cannam@95 201 if (rnk == 1 && X(is_prime)(sz->dims[0].n))
cannam@95 202 sz->dims[0].b[IB] = sz->dims[0].b[OB] = sz->dims[0].n;
cannam@95 203
cannam@95 204 XM(dtensor_destroy)(sz0);
cannam@95 205 sz0 = XM(dtensor_canonical)(sz, 0);
cannam@95 206 XM(dtensor_destroy)(sz);
cannam@95 207 return sz0;
cannam@95 208 }
cannam@95 209
cannam@95 210 /* allocate simple local (serial) dims array corresponding to n[rnk] */
cannam@95 211 static XM(ddim) *simple_dims(int rnk, const ptrdiff_t *n)
cannam@95 212 {
cannam@95 213 XM(ddim) *dims = (XM(ddim) *) MALLOC(sizeof(XM(ddim)) * rnk,
cannam@95 214 TENSORS);
cannam@95 215 int i;
cannam@95 216 for (i = 0; i < rnk; ++i)
cannam@95 217 dims[i].n = dims[i].ib = dims[i].ob = n[i];
cannam@95 218 return dims;
cannam@95 219 }
cannam@95 220
cannam@95 221 /*************************************************************************/
cannam@95 222
cannam@95 223 static void local_size(int my_pe, const dtensor *sz, block_kind k,
cannam@95 224 ptrdiff_t *local_n, ptrdiff_t *local_start)
cannam@95 225 {
cannam@95 226 int i;
cannam@95 227 if (my_pe >= XM(num_blocks_total)(sz, k))
cannam@95 228 for (i = 0; i < sz->rnk; ++i)
cannam@95 229 local_n[i] = local_start[i] = 0;
cannam@95 230 else {
cannam@95 231 XM(block_coords)(sz, k, my_pe, local_start);
cannam@95 232 for (i = 0; i < sz->rnk; ++i) {
cannam@95 233 local_n[i] = XM(block)(sz->dims[i].n, sz->dims[i].b[k],
cannam@95 234 local_start[i]);
cannam@95 235 local_start[i] *= sz->dims[i].b[k];
cannam@95 236 }
cannam@95 237 }
cannam@95 238 }
cannam@95 239
cannam@95 240 static INT prod(int rnk, const ptrdiff_t *local_n)
cannam@95 241 {
cannam@95 242 int i;
cannam@95 243 INT N = 1;
cannam@95 244 for (i = 0; i < rnk; ++i) N *= local_n[i];
cannam@95 245 return N;
cannam@95 246 }
cannam@95 247
cannam@95 248 ptrdiff_t XM(local_size_guru)(int rnk, const XM(ddim) *dims0,
cannam@95 249 ptrdiff_t howmany, MPI_Comm comm,
cannam@95 250 ptrdiff_t *local_n_in,
cannam@95 251 ptrdiff_t *local_start_in,
cannam@95 252 ptrdiff_t *local_n_out,
cannam@95 253 ptrdiff_t *local_start_out,
cannam@95 254 int sign, unsigned flags)
cannam@95 255 {
cannam@95 256 INT N;
cannam@95 257 int my_pe, n_pes, i;
cannam@95 258 dtensor *sz;
cannam@95 259
cannam@95 260 if (rnk == 0)
cannam@95 261 return howmany;
cannam@95 262
cannam@95 263 MPI_Comm_rank(comm, &my_pe);
cannam@95 264 MPI_Comm_size(comm, &n_pes);
cannam@95 265 sz = default_sz(rnk, dims0, n_pes, 0);
cannam@95 266
cannam@95 267 /* Now, we must figure out how much local space the user should
cannam@95 268 allocate (or at least an upper bound). This depends strongly
cannam@95 269 on the exact algorithms we employ...ugh! FIXME: get this info
cannam@95 270 from the solvers somehow? */
cannam@95 271 N = 1; /* never return zero allocation size */
cannam@95 272 if (rnk > 1 && XM(is_block1d)(sz, IB) && XM(is_block1d)(sz, OB)) {
cannam@95 273 INT Nafter;
cannam@95 274 ddim odims[2];
cannam@95 275
cannam@95 276 /* dft-rank-geq2-transposed */
cannam@95 277 odims[0] = sz->dims[0]; odims[1] = sz->dims[1]; /* save */
cannam@95 278 /* we may need extra space for transposed intermediate data */
cannam@95 279 for (i = 0; i < 2; ++i)
cannam@95 280 if (XM(num_blocks)(sz->dims[i].n, sz->dims[i].b[IB]) == 1 &&
cannam@95 281 XM(num_blocks)(sz->dims[i].n, sz->dims[i].b[OB]) == 1) {
cannam@95 282 sz->dims[i].b[IB]
cannam@95 283 = XM(default_block)(sz->dims[i].n, n_pes);
cannam@95 284 sz->dims[1-i].b[IB] = sz->dims[1-i].n;
cannam@95 285 local_size(my_pe, sz, IB, local_n_in, local_start_in);
cannam@95 286 N = X(imax)(N, prod(rnk, local_n_in));
cannam@95 287 sz->dims[i] = odims[i];
cannam@95 288 sz->dims[1-i] = odims[1-i];
cannam@95 289 break;
cannam@95 290 }
cannam@95 291
cannam@95 292 /* dft-rank-geq2 */
cannam@95 293 Nafter = howmany;
cannam@95 294 for (i = 1; i < sz->rnk; ++i) Nafter *= sz->dims[i].n;
cannam@95 295 N = X(imax)(N, (sz->dims[0].n
cannam@95 296 * XM(block)(Nafter, XM(default_block)(Nafter, n_pes),
cannam@95 297 my_pe) + howmany - 1) / howmany);
cannam@95 298
cannam@95 299 /* dft-rank-geq2 with dimensions swapped */
cannam@95 300 Nafter = howmany * sz->dims[0].n;
cannam@95 301 for (i = 2; i < sz->rnk; ++i) Nafter *= sz->dims[i].n;
cannam@95 302 N = X(imax)(N, (sz->dims[1].n
cannam@95 303 * XM(block)(Nafter, XM(default_block)(Nafter, n_pes),
cannam@95 304 my_pe) + howmany - 1) / howmany);
cannam@95 305 }
cannam@95 306 else if (rnk == 1) {
cannam@95 307 if (howmany >= n_pes && !MPI_FLAGS(flags)) { /* dft-rank1-bigvec */
cannam@95 308 ptrdiff_t n[2], start[2];
cannam@95 309 dtensor *sz2 = XM(mkdtensor)(2);
cannam@95 310 sz2->dims[0] = sz->dims[0];
cannam@95 311 sz2->dims[0].b[IB] = sz->dims[0].n;
cannam@95 312 sz2->dims[1].n = sz2->dims[1].b[OB] = howmany;
cannam@95 313 sz2->dims[1].b[IB] = XM(default_block)(howmany, n_pes);
cannam@95 314 local_size(my_pe, sz2, IB, n, start);
cannam@95 315 XM(dtensor_destroy)(sz2);
cannam@95 316 N = X(imax)(N, (prod(2, n) + howmany - 1) / howmany);
cannam@95 317 }
cannam@95 318 else { /* dft-rank1 */
cannam@95 319 INT r, m, rblock[2], mblock[2];
cannam@95 320
cannam@95 321 /* Since the 1d transforms are so different, we require
cannam@95 322 the user to call local_size_1d for this case. Ugh. */
cannam@95 323 CK(sign == FFTW_FORWARD || sign == FFTW_BACKWARD);
cannam@95 324
cannam@95 325 if ((r = XM(choose_radix)(sz->dims[0], n_pes, flags, sign,
cannam@95 326 rblock, mblock))) {
cannam@95 327 m = sz->dims[0].n / r;
cannam@95 328 if (flags & FFTW_MPI_SCRAMBLED_IN)
cannam@95 329 sz->dims[0].b[IB] = rblock[IB] * m;
cannam@95 330 else { /* !SCRAMBLED_IN */
cannam@95 331 sz->dims[0].b[IB] = r * mblock[IB];
cannam@95 332 N = X(imax)(N, rblock[IB] * m);
cannam@95 333 }
cannam@95 334 if (flags & FFTW_MPI_SCRAMBLED_OUT)
cannam@95 335 sz->dims[0].b[OB] = r * mblock[OB];
cannam@95 336 else { /* !SCRAMBLED_OUT */
cannam@95 337 N = X(imax)(N, r * mblock[OB]);
cannam@95 338 sz->dims[0].b[OB] = rblock[OB] * m;
cannam@95 339 }
cannam@95 340 }
cannam@95 341 }
cannam@95 342 }
cannam@95 343
cannam@95 344 local_size(my_pe, sz, IB, local_n_in, local_start_in);
cannam@95 345 local_size(my_pe, sz, OB, local_n_out, local_start_out);
cannam@95 346
cannam@95 347 /* at least, make sure we have enough space to store input & output */
cannam@95 348 N = X(imax)(N, X(imax)(prod(rnk, local_n_in), prod(rnk, local_n_out)));
cannam@95 349
cannam@95 350 XM(dtensor_destroy)(sz);
cannam@95 351 return N * howmany;
cannam@95 352 }
cannam@95 353
cannam@95 354 ptrdiff_t XM(local_size_many_transposed)(int rnk, const ptrdiff_t *n,
cannam@95 355 ptrdiff_t howmany,
cannam@95 356 ptrdiff_t xblock, ptrdiff_t yblock,
cannam@95 357 MPI_Comm comm,
cannam@95 358 ptrdiff_t *local_nx,
cannam@95 359 ptrdiff_t *local_x_start,
cannam@95 360 ptrdiff_t *local_ny,
cannam@95 361 ptrdiff_t *local_y_start)
cannam@95 362 {
cannam@95 363 ptrdiff_t N;
cannam@95 364 XM(ddim) *dims;
cannam@95 365 ptrdiff_t *local;
cannam@95 366
cannam@95 367 if (rnk == 0) {
cannam@95 368 *local_nx = *local_ny = 1;
cannam@95 369 *local_x_start = *local_y_start = 0;
cannam@95 370 return howmany;
cannam@95 371 }
cannam@95 372
cannam@95 373 dims = simple_dims(rnk, n);
cannam@95 374 local = (ptrdiff_t *) MALLOC(sizeof(ptrdiff_t) * rnk * 4, TENSORS);
cannam@95 375
cannam@95 376 /* default 1d block distribution, with transposed output
cannam@95 377 if yblock < n[1] */
cannam@95 378 dims[0].ib = xblock;
cannam@95 379 if (rnk > 1) {
cannam@95 380 if (yblock < n[1])
cannam@95 381 dims[1].ob = yblock;
cannam@95 382 else
cannam@95 383 dims[0].ob = xblock;
cannam@95 384 }
cannam@95 385 else
cannam@95 386 dims[0].ob = xblock; /* FIXME: 1d not really supported here
cannam@95 387 since we don't have flags/sign */
cannam@95 388
cannam@95 389 N = XM(local_size_guru)(rnk, dims, howmany, comm,
cannam@95 390 local, local + rnk,
cannam@95 391 local + 2*rnk, local + 3*rnk,
cannam@95 392 0, 0);
cannam@95 393 *local_nx = local[0];
cannam@95 394 *local_x_start = local[rnk];
cannam@95 395 if (rnk > 1) {
cannam@95 396 *local_ny = local[2*rnk + 1];
cannam@95 397 *local_y_start = local[3*rnk + 1];
cannam@95 398 }
cannam@95 399 else {
cannam@95 400 *local_ny = *local_nx;
cannam@95 401 *local_y_start = *local_x_start;
cannam@95 402 }
cannam@95 403 X(ifree)(local);
cannam@95 404 X(ifree)(dims);
cannam@95 405 return N;
cannam@95 406 }
cannam@95 407
cannam@95 408 ptrdiff_t XM(local_size_many)(int rnk, const ptrdiff_t *n,
cannam@95 409 ptrdiff_t howmany,
cannam@95 410 ptrdiff_t xblock,
cannam@95 411 MPI_Comm comm,
cannam@95 412 ptrdiff_t *local_nx,
cannam@95 413 ptrdiff_t *local_x_start)
cannam@95 414 {
cannam@95 415 ptrdiff_t local_ny, local_y_start;
cannam@95 416 return XM(local_size_many_transposed)(rnk, n, howmany,
cannam@95 417 xblock, rnk > 1
cannam@95 418 ? n[1] : FFTW_MPI_DEFAULT_BLOCK,
cannam@95 419 comm,
cannam@95 420 local_nx, local_x_start,
cannam@95 421 &local_ny, &local_y_start);
cannam@95 422 }
cannam@95 423
cannam@95 424
cannam@95 425 ptrdiff_t XM(local_size_transposed)(int rnk, const ptrdiff_t *n,
cannam@95 426 MPI_Comm comm,
cannam@95 427 ptrdiff_t *local_nx,
cannam@95 428 ptrdiff_t *local_x_start,
cannam@95 429 ptrdiff_t *local_ny,
cannam@95 430 ptrdiff_t *local_y_start)
cannam@95 431 {
cannam@95 432 return XM(local_size_many_transposed)(rnk, n, 1,
cannam@95 433 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 434 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 435 comm,
cannam@95 436 local_nx, local_x_start,
cannam@95 437 local_ny, local_y_start);
cannam@95 438 }
cannam@95 439
cannam@95 440 ptrdiff_t XM(local_size)(int rnk, const ptrdiff_t *n,
cannam@95 441 MPI_Comm comm,
cannam@95 442 ptrdiff_t *local_nx,
cannam@95 443 ptrdiff_t *local_x_start)
cannam@95 444 {
cannam@95 445 return XM(local_size_many)(rnk, n, 1, FFTW_MPI_DEFAULT_BLOCK, comm,
cannam@95 446 local_nx, local_x_start);
cannam@95 447 }
cannam@95 448
cannam@95 449 ptrdiff_t XM(local_size_many_1d)(ptrdiff_t nx, ptrdiff_t howmany,
cannam@95 450 MPI_Comm comm, int sign, unsigned flags,
cannam@95 451 ptrdiff_t *local_nx, ptrdiff_t *local_x_start,
cannam@95 452 ptrdiff_t *local_ny, ptrdiff_t *local_y_start)
cannam@95 453 {
cannam@95 454 XM(ddim) d;
cannam@95 455 d.n = nx;
cannam@95 456 d.ib = d.ob = FFTW_MPI_DEFAULT_BLOCK;
cannam@95 457 return XM(local_size_guru)(1, &d, howmany, comm,
cannam@95 458 local_nx, local_x_start,
cannam@95 459 local_ny, local_y_start, sign, flags);
cannam@95 460 }
cannam@95 461
cannam@95 462 ptrdiff_t XM(local_size_1d)(ptrdiff_t nx,
cannam@95 463 MPI_Comm comm, int sign, unsigned flags,
cannam@95 464 ptrdiff_t *local_nx, ptrdiff_t *local_x_start,
cannam@95 465 ptrdiff_t *local_ny, ptrdiff_t *local_y_start)
cannam@95 466 {
cannam@95 467 return XM(local_size_many_1d)(nx, 1, comm, sign, flags,
cannam@95 468 local_nx, local_x_start,
cannam@95 469 local_ny, local_y_start);
cannam@95 470 }
cannam@95 471
cannam@95 472 ptrdiff_t XM(local_size_2d_transposed)(ptrdiff_t nx, ptrdiff_t ny,
cannam@95 473 MPI_Comm comm,
cannam@95 474 ptrdiff_t *local_nx,
cannam@95 475 ptrdiff_t *local_x_start,
cannam@95 476 ptrdiff_t *local_ny,
cannam@95 477 ptrdiff_t *local_y_start)
cannam@95 478 {
cannam@95 479 ptrdiff_t n[2];
cannam@95 480 n[0] = nx; n[1] = ny;
cannam@95 481 return XM(local_size_transposed)(2, n, comm,
cannam@95 482 local_nx, local_x_start,
cannam@95 483 local_ny, local_y_start);
cannam@95 484 }
cannam@95 485
cannam@95 486 ptrdiff_t XM(local_size_2d)(ptrdiff_t nx, ptrdiff_t ny, MPI_Comm comm,
cannam@95 487 ptrdiff_t *local_nx, ptrdiff_t *local_x_start)
cannam@95 488 {
cannam@95 489 ptrdiff_t n[2];
cannam@95 490 n[0] = nx; n[1] = ny;
cannam@95 491 return XM(local_size)(2, n, comm, local_nx, local_x_start);
cannam@95 492 }
cannam@95 493
cannam@95 494 ptrdiff_t XM(local_size_3d_transposed)(ptrdiff_t nx, ptrdiff_t ny,
cannam@95 495 ptrdiff_t nz,
cannam@95 496 MPI_Comm comm,
cannam@95 497 ptrdiff_t *local_nx,
cannam@95 498 ptrdiff_t *local_x_start,
cannam@95 499 ptrdiff_t *local_ny,
cannam@95 500 ptrdiff_t *local_y_start)
cannam@95 501 {
cannam@95 502 ptrdiff_t n[3];
cannam@95 503 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 504 return XM(local_size_transposed)(3, n, comm,
cannam@95 505 local_nx, local_x_start,
cannam@95 506 local_ny, local_y_start);
cannam@95 507 }
cannam@95 508
cannam@95 509 ptrdiff_t XM(local_size_3d)(ptrdiff_t nx, ptrdiff_t ny, ptrdiff_t nz,
cannam@95 510 MPI_Comm comm,
cannam@95 511 ptrdiff_t *local_nx, ptrdiff_t *local_x_start)
cannam@95 512 {
cannam@95 513 ptrdiff_t n[3];
cannam@95 514 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 515 return XM(local_size)(3, n, comm, local_nx, local_x_start);
cannam@95 516 }
cannam@95 517
cannam@95 518 /*************************************************************************/
cannam@95 519 /* Transpose API */
cannam@95 520
cannam@95 521 X(plan) XM(plan_many_transpose)(ptrdiff_t nx, ptrdiff_t ny,
cannam@95 522 ptrdiff_t howmany,
cannam@95 523 ptrdiff_t xblock, ptrdiff_t yblock,
cannam@95 524 R *in, R *out,
cannam@95 525 MPI_Comm comm, unsigned flags)
cannam@95 526 {
cannam@95 527 int n_pes;
cannam@95 528 XM(init)();
cannam@95 529
cannam@95 530 if (howmany < 0 || xblock < 0 || yblock < 0 ||
cannam@95 531 nx <= 0 || ny <= 0) return 0;
cannam@95 532
cannam@95 533 MPI_Comm_size(comm, &n_pes);
cannam@95 534 if (!xblock) xblock = XM(default_block)(nx, n_pes);
cannam@95 535 if (!yblock) yblock = XM(default_block)(ny, n_pes);
cannam@95 536 if (n_pes < XM(num_blocks)(nx, xblock)
cannam@95 537 || n_pes < XM(num_blocks)(ny, yblock))
cannam@95 538 return 0;
cannam@95 539
cannam@95 540 return
cannam@95 541 X(mkapiplan)(FFTW_FORWARD, flags,
cannam@95 542 XM(mkproblem_transpose)(nx, ny, howmany,
cannam@95 543 in, out, xblock, yblock,
cannam@95 544 comm, MPI_FLAGS(flags)));
cannam@95 545 }
cannam@95 546
cannam@95 547 X(plan) XM(plan_transpose)(ptrdiff_t nx, ptrdiff_t ny, R *in, R *out,
cannam@95 548 MPI_Comm comm, unsigned flags)
cannam@95 549
cannam@95 550 {
cannam@95 551 return XM(plan_many_transpose)(nx, ny, 1,
cannam@95 552 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 553 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 554 in, out, comm, flags);
cannam@95 555 }
cannam@95 556
cannam@95 557 /*************************************************************************/
cannam@95 558 /* Complex DFT API */
cannam@95 559
cannam@95 560 X(plan) XM(plan_guru_dft)(int rnk, const XM(ddim) *dims0,
cannam@95 561 ptrdiff_t howmany,
cannam@95 562 C *in, C *out,
cannam@95 563 MPI_Comm comm, int sign, unsigned flags)
cannam@95 564 {
cannam@95 565 int n_pes, i;
cannam@95 566 dtensor *sz;
cannam@95 567
cannam@95 568 XM(init)();
cannam@95 569
cannam@95 570 if (howmany < 0 || rnk < 1) return 0;
cannam@95 571 for (i = 0; i < rnk; ++i)
cannam@95 572 if (dims0[i].n < 1 || dims0[i].ib < 0 || dims0[i].ob < 0)
cannam@95 573 return 0;
cannam@95 574
cannam@95 575 MPI_Comm_size(comm, &n_pes);
cannam@95 576 sz = default_sz(rnk, dims0, n_pes, 0);
cannam@95 577
cannam@95 578 if (XM(num_blocks_total)(sz, IB) > n_pes
cannam@95 579 || XM(num_blocks_total)(sz, OB) > n_pes) {
cannam@95 580 XM(dtensor_destroy)(sz);
cannam@95 581 return 0;
cannam@95 582 }
cannam@95 583
cannam@95 584 return
cannam@95 585 X(mkapiplan)(sign, flags,
cannam@95 586 XM(mkproblem_dft_d)(sz, howmany,
cannam@95 587 (R *) in, (R *) out,
cannam@95 588 comm, sign,
cannam@95 589 MPI_FLAGS(flags)));
cannam@95 590 }
cannam@95 591
cannam@95 592 X(plan) XM(plan_many_dft)(int rnk, const ptrdiff_t *n,
cannam@95 593 ptrdiff_t howmany,
cannam@95 594 ptrdiff_t iblock, ptrdiff_t oblock,
cannam@95 595 C *in, C *out,
cannam@95 596 MPI_Comm comm, int sign, unsigned flags)
cannam@95 597 {
cannam@95 598 XM(ddim) *dims = simple_dims(rnk, n);
cannam@95 599 X(plan) pln;
cannam@95 600
cannam@95 601 if (rnk == 1) {
cannam@95 602 dims[0].ib = iblock;
cannam@95 603 dims[0].ob = oblock;
cannam@95 604 }
cannam@95 605 else if (rnk > 1) {
cannam@95 606 dims[0 != (flags & FFTW_MPI_TRANSPOSED_IN)].ib = iblock;
cannam@95 607 dims[0 != (flags & FFTW_MPI_TRANSPOSED_OUT)].ob = oblock;
cannam@95 608 }
cannam@95 609
cannam@95 610 pln = XM(plan_guru_dft)(rnk,dims,howmany, in,out, comm, sign, flags);
cannam@95 611 X(ifree)(dims);
cannam@95 612 return pln;
cannam@95 613 }
cannam@95 614
cannam@95 615 X(plan) XM(plan_dft)(int rnk, const ptrdiff_t *n, C *in, C *out,
cannam@95 616 MPI_Comm comm, int sign, unsigned flags)
cannam@95 617 {
cannam@95 618 return XM(plan_many_dft)(rnk, n, 1,
cannam@95 619 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 620 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 621 in, out, comm, sign, flags);
cannam@95 622 }
cannam@95 623
cannam@95 624 X(plan) XM(plan_dft_1d)(ptrdiff_t nx, C *in, C *out,
cannam@95 625 MPI_Comm comm, int sign, unsigned flags)
cannam@95 626 {
cannam@95 627 return XM(plan_dft)(1, &nx, in, out, comm, sign, flags);
cannam@95 628 }
cannam@95 629
cannam@95 630 X(plan) XM(plan_dft_2d)(ptrdiff_t nx, ptrdiff_t ny, C *in, C *out,
cannam@95 631 MPI_Comm comm, int sign, unsigned flags)
cannam@95 632 {
cannam@95 633 ptrdiff_t n[2];
cannam@95 634 n[0] = nx; n[1] = ny;
cannam@95 635 return XM(plan_dft)(2, n, in, out, comm, sign, flags);
cannam@95 636 }
cannam@95 637
cannam@95 638 X(plan) XM(plan_dft_3d)(ptrdiff_t nx, ptrdiff_t ny, ptrdiff_t nz,
cannam@95 639 C *in, C *out,
cannam@95 640 MPI_Comm comm, int sign, unsigned flags)
cannam@95 641 {
cannam@95 642 ptrdiff_t n[3];
cannam@95 643 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 644 return XM(plan_dft)(3, n, in, out, comm, sign, flags);
cannam@95 645 }
cannam@95 646
cannam@95 647 /*************************************************************************/
cannam@95 648 /* R2R API */
cannam@95 649
cannam@95 650 X(plan) XM(plan_guru_r2r)(int rnk, const XM(ddim) *dims0,
cannam@95 651 ptrdiff_t howmany,
cannam@95 652 R *in, R *out,
cannam@95 653 MPI_Comm comm, const X(r2r_kind) *kind,
cannam@95 654 unsigned flags)
cannam@95 655 {
cannam@95 656 int n_pes, i;
cannam@95 657 dtensor *sz;
cannam@95 658 rdft_kind *k;
cannam@95 659 X(plan) pln;
cannam@95 660
cannam@95 661 XM(init)();
cannam@95 662
cannam@95 663 if (howmany < 0 || rnk < 1) return 0;
cannam@95 664 for (i = 0; i < rnk; ++i)
cannam@95 665 if (dims0[i].n < 1 || dims0[i].ib < 0 || dims0[i].ob < 0)
cannam@95 666 return 0;
cannam@95 667
cannam@95 668 k = X(map_r2r_kind)(rnk, kind);
cannam@95 669
cannam@95 670 MPI_Comm_size(comm, &n_pes);
cannam@95 671 sz = default_sz(rnk, dims0, n_pes, 0);
cannam@95 672
cannam@95 673 if (XM(num_blocks_total)(sz, IB) > n_pes
cannam@95 674 || XM(num_blocks_total)(sz, OB) > n_pes) {
cannam@95 675 XM(dtensor_destroy)(sz);
cannam@95 676 return 0;
cannam@95 677 }
cannam@95 678
cannam@95 679 pln = X(mkapiplan)(0, flags,
cannam@95 680 XM(mkproblem_rdft_d)(sz, howmany,
cannam@95 681 in, out,
cannam@95 682 comm, k, MPI_FLAGS(flags)));
cannam@95 683 X(ifree0)(k);
cannam@95 684 return pln;
cannam@95 685 }
cannam@95 686
cannam@95 687 X(plan) XM(plan_many_r2r)(int rnk, const ptrdiff_t *n,
cannam@95 688 ptrdiff_t howmany,
cannam@95 689 ptrdiff_t iblock, ptrdiff_t oblock,
cannam@95 690 R *in, R *out,
cannam@95 691 MPI_Comm comm, const X(r2r_kind) *kind,
cannam@95 692 unsigned flags)
cannam@95 693 {
cannam@95 694 XM(ddim) *dims = simple_dims(rnk, n);
cannam@95 695 X(plan) pln;
cannam@95 696
cannam@95 697 if (rnk == 1) {
cannam@95 698 dims[0].ib = iblock;
cannam@95 699 dims[0].ob = oblock;
cannam@95 700 }
cannam@95 701 else if (rnk > 1) {
cannam@95 702 dims[0 != (flags & FFTW_MPI_TRANSPOSED_IN)].ib = iblock;
cannam@95 703 dims[0 != (flags & FFTW_MPI_TRANSPOSED_OUT)].ob = oblock;
cannam@95 704 }
cannam@95 705
cannam@95 706 pln = XM(plan_guru_r2r)(rnk,dims,howmany, in,out, comm, kind, flags);
cannam@95 707 X(ifree)(dims);
cannam@95 708 return pln;
cannam@95 709 }
cannam@95 710
cannam@95 711 X(plan) XM(plan_r2r)(int rnk, const ptrdiff_t *n, R *in, R *out,
cannam@95 712 MPI_Comm comm,
cannam@95 713 const X(r2r_kind) *kind,
cannam@95 714 unsigned flags)
cannam@95 715 {
cannam@95 716 return XM(plan_many_r2r)(rnk, n, 1,
cannam@95 717 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 718 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 719 in, out, comm, kind, flags);
cannam@95 720 }
cannam@95 721
cannam@95 722 X(plan) XM(plan_r2r_2d)(ptrdiff_t nx, ptrdiff_t ny, R *in, R *out,
cannam@95 723 MPI_Comm comm,
cannam@95 724 X(r2r_kind) kindx, X(r2r_kind) kindy,
cannam@95 725 unsigned flags)
cannam@95 726 {
cannam@95 727 ptrdiff_t n[2];
cannam@95 728 X(r2r_kind) kind[2];
cannam@95 729 n[0] = nx; n[1] = ny;
cannam@95 730 kind[0] = kindx; kind[1] = kindy;
cannam@95 731 return XM(plan_r2r)(2, n, in, out, comm, kind, flags);
cannam@95 732 }
cannam@95 733
cannam@95 734 X(plan) XM(plan_r2r_3d)(ptrdiff_t nx, ptrdiff_t ny, ptrdiff_t nz,
cannam@95 735 R *in, R *out,
cannam@95 736 MPI_Comm comm,
cannam@95 737 X(r2r_kind) kindx, X(r2r_kind) kindy,
cannam@95 738 X(r2r_kind) kindz,
cannam@95 739 unsigned flags)
cannam@95 740 {
cannam@95 741 ptrdiff_t n[3];
cannam@95 742 X(r2r_kind) kind[3];
cannam@95 743 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 744 kind[0] = kindx; kind[1] = kindy; kind[2] = kindz;
cannam@95 745 return XM(plan_r2r)(3, n, in, out, comm, kind, flags);
cannam@95 746 }
cannam@95 747
cannam@95 748 /*************************************************************************/
cannam@95 749 /* R2C/C2R API */
cannam@95 750
cannam@95 751 static X(plan) plan_guru_rdft2(int rnk, const XM(ddim) *dims0,
cannam@95 752 ptrdiff_t howmany,
cannam@95 753 R *r, C *c,
cannam@95 754 MPI_Comm comm, rdft_kind kind, unsigned flags)
cannam@95 755 {
cannam@95 756 int n_pes, i;
cannam@95 757 dtensor *sz;
cannam@95 758 R *cr = (R *) c;
cannam@95 759
cannam@95 760 XM(init)();
cannam@95 761
cannam@95 762 if (howmany < 0 || rnk < 2) return 0;
cannam@95 763 for (i = 0; i < rnk; ++i)
cannam@95 764 if (dims0[i].n < 1 || dims0[i].ib < 0 || dims0[i].ob < 0)
cannam@95 765 return 0;
cannam@95 766
cannam@95 767 MPI_Comm_size(comm, &n_pes);
cannam@95 768 sz = default_sz(rnk, dims0, n_pes, 1);
cannam@95 769
cannam@95 770 sz->dims[rnk-1].n = dims0[rnk-1].n / 2 + 1;
cannam@95 771 if (XM(num_blocks_total)(sz, IB) > n_pes
cannam@95 772 || XM(num_blocks_total)(sz, OB) > n_pes) {
cannam@95 773 XM(dtensor_destroy)(sz);
cannam@95 774 return 0;
cannam@95 775 }
cannam@95 776 sz->dims[rnk-1].n = dims0[rnk-1].n;
cannam@95 777
cannam@95 778 if (kind == R2HC)
cannam@95 779 return X(mkapiplan)(0, flags,
cannam@95 780 XM(mkproblem_rdft2_d)(sz, howmany,
cannam@95 781 r, cr, comm, R2HC,
cannam@95 782 MPI_FLAGS(flags)));
cannam@95 783 else
cannam@95 784 return X(mkapiplan)(0, flags,
cannam@95 785 XM(mkproblem_rdft2_d)(sz, howmany,
cannam@95 786 cr, r, comm, HC2R,
cannam@95 787 MPI_FLAGS(flags)));
cannam@95 788 }
cannam@95 789
cannam@95 790 X(plan) XM(plan_many_dft_r2c)(int rnk, const ptrdiff_t *n,
cannam@95 791 ptrdiff_t howmany,
cannam@95 792 ptrdiff_t iblock, ptrdiff_t oblock,
cannam@95 793 R *in, C *out,
cannam@95 794 MPI_Comm comm, unsigned flags)
cannam@95 795 {
cannam@95 796 XM(ddim) *dims = simple_dims(rnk, n);
cannam@95 797 X(plan) pln;
cannam@95 798
cannam@95 799 if (rnk == 1) {
cannam@95 800 dims[0].ib = iblock;
cannam@95 801 dims[0].ob = oblock;
cannam@95 802 }
cannam@95 803 else if (rnk > 1) {
cannam@95 804 dims[0 != (flags & FFTW_MPI_TRANSPOSED_IN)].ib = iblock;
cannam@95 805 dims[0 != (flags & FFTW_MPI_TRANSPOSED_OUT)].ob = oblock;
cannam@95 806 }
cannam@95 807
cannam@95 808 pln = plan_guru_rdft2(rnk,dims,howmany, in,out, comm, R2HC, flags);
cannam@95 809 X(ifree)(dims);
cannam@95 810 return pln;
cannam@95 811 }
cannam@95 812
cannam@95 813 X(plan) XM(plan_many_dft_c2r)(int rnk, const ptrdiff_t *n,
cannam@95 814 ptrdiff_t howmany,
cannam@95 815 ptrdiff_t iblock, ptrdiff_t oblock,
cannam@95 816 C *in, R *out,
cannam@95 817 MPI_Comm comm, unsigned flags)
cannam@95 818 {
cannam@95 819 XM(ddim) *dims = simple_dims(rnk, n);
cannam@95 820 X(plan) pln;
cannam@95 821
cannam@95 822 if (rnk == 1) {
cannam@95 823 dims[0].ib = iblock;
cannam@95 824 dims[0].ob = oblock;
cannam@95 825 }
cannam@95 826 else if (rnk > 1) {
cannam@95 827 dims[0 != (flags & FFTW_MPI_TRANSPOSED_IN)].ib = iblock;
cannam@95 828 dims[0 != (flags & FFTW_MPI_TRANSPOSED_OUT)].ob = oblock;
cannam@95 829 }
cannam@95 830
cannam@95 831 pln = plan_guru_rdft2(rnk,dims,howmany, out,in, comm, HC2R, flags);
cannam@95 832 X(ifree)(dims);
cannam@95 833 return pln;
cannam@95 834 }
cannam@95 835
cannam@95 836 X(plan) XM(plan_dft_r2c)(int rnk, const ptrdiff_t *n, R *in, C *out,
cannam@95 837 MPI_Comm comm, unsigned flags)
cannam@95 838 {
cannam@95 839 return XM(plan_many_dft_r2c)(rnk, n, 1,
cannam@95 840 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 841 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 842 in, out, comm, flags);
cannam@95 843 }
cannam@95 844
cannam@95 845 X(plan) XM(plan_dft_r2c_2d)(ptrdiff_t nx, ptrdiff_t ny, R *in, C *out,
cannam@95 846 MPI_Comm comm, unsigned flags)
cannam@95 847 {
cannam@95 848 ptrdiff_t n[2];
cannam@95 849 n[0] = nx; n[1] = ny;
cannam@95 850 return XM(plan_dft_r2c)(2, n, in, out, comm, flags);
cannam@95 851 }
cannam@95 852
cannam@95 853 X(plan) XM(plan_dft_r2c_3d)(ptrdiff_t nx, ptrdiff_t ny, ptrdiff_t nz,
cannam@95 854 R *in, C *out, MPI_Comm comm, unsigned flags)
cannam@95 855 {
cannam@95 856 ptrdiff_t n[3];
cannam@95 857 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 858 return XM(plan_dft_r2c)(3, n, in, out, comm, flags);
cannam@95 859 }
cannam@95 860
cannam@95 861 X(plan) XM(plan_dft_c2r)(int rnk, const ptrdiff_t *n, C *in, R *out,
cannam@95 862 MPI_Comm comm, unsigned flags)
cannam@95 863 {
cannam@95 864 return XM(plan_many_dft_c2r)(rnk, n, 1,
cannam@95 865 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 866 FFTW_MPI_DEFAULT_BLOCK,
cannam@95 867 in, out, comm, flags);
cannam@95 868 }
cannam@95 869
cannam@95 870 X(plan) XM(plan_dft_c2r_2d)(ptrdiff_t nx, ptrdiff_t ny, C *in, R *out,
cannam@95 871 MPI_Comm comm, unsigned flags)
cannam@95 872 {
cannam@95 873 ptrdiff_t n[2];
cannam@95 874 n[0] = nx; n[1] = ny;
cannam@95 875 return XM(plan_dft_c2r)(2, n, in, out, comm, flags);
cannam@95 876 }
cannam@95 877
cannam@95 878 X(plan) XM(plan_dft_c2r_3d)(ptrdiff_t nx, ptrdiff_t ny, ptrdiff_t nz,
cannam@95 879 C *in, R *out, MPI_Comm comm, unsigned flags)
cannam@95 880 {
cannam@95 881 ptrdiff_t n[3];
cannam@95 882 n[0] = nx; n[1] = ny; n[2] = nz;
cannam@95 883 return XM(plan_dft_c2r)(3, n, in, out, comm, flags);
cannam@95 884 }
cannam@95 885
cannam@95 886 /*************************************************************************/
cannam@95 887 /* New-array execute functions */
cannam@95 888
cannam@95 889 void XM(execute_dft)(const X(plan) p, C *in, C *out) {
cannam@95 890 /* internally, MPI plans are just rdft plans */
cannam@95 891 X(execute_r2r)(p, (R*) in, (R*) out);
cannam@95 892 }
cannam@95 893
cannam@95 894 void XM(execute_dft_r2c)(const X(plan) p, R *in, C *out) {
cannam@95 895 /* internally, MPI plans are just rdft plans */
cannam@95 896 X(execute_r2r)(p, in, (R*) out);
cannam@95 897 }
cannam@95 898
cannam@95 899 void XM(execute_dft_c2r)(const X(plan) p, C *in, R *out) {
cannam@95 900 /* internally, MPI plans are just rdft plans */
cannam@95 901 X(execute_r2r)(p, (R*) in, out);
cannam@95 902 }
cannam@95 903
cannam@95 904 void XM(execute_r2r)(const X(plan) p, R *in, R *out) {
cannam@95 905 /* internally, MPI plans are just rdft plans */
cannam@95 906 X(execute_r2r)(p, in, out);
cannam@95 907 }