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