annotate src/fftw-3.3.8/mpi/transpose-pairwise.c @ 82:d0c2a83c1364

Add FFTW 3.3.8 source, and a Linux build
author Chris Cannam
date Tue, 19 Nov 2019 14:52:55 +0000
parents
children
rev   line source
Chris@82 1 /*
Chris@82 2 * Copyright (c) 2003, 2007-14 Matteo Frigo
Chris@82 3 * Copyright (c) 2003, 2007-14 Massachusetts Institute of Technology
Chris@82 4 *
Chris@82 5 * This program is free software; you can redistribute it and/or modify
Chris@82 6 * it under the terms of the GNU General Public License as published by
Chris@82 7 * the Free Software Foundation; either version 2 of the License, or
Chris@82 8 * (at your option) any later version.
Chris@82 9 *
Chris@82 10 * This program is distributed in the hope that it will be useful,
Chris@82 11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
Chris@82 12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
Chris@82 13 * GNU General Public License for more details.
Chris@82 14 *
Chris@82 15 * You should have received a copy of the GNU General Public License
Chris@82 16 * along with this program; if not, write to the Free Software
Chris@82 17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
Chris@82 18 *
Chris@82 19 */
Chris@82 20
Chris@82 21 /* Distributed transposes using a sequence of carefully scheduled
Chris@82 22 pairwise exchanges. This has the advantage that it can be done
Chris@82 23 in-place, or out-of-place while preserving the input, using buffer
Chris@82 24 space proportional to the local size divided by the number of
Chris@82 25 processes (i.e. to the total array size divided by the number of
Chris@82 26 processes squared). */
Chris@82 27
Chris@82 28 #include "mpi-transpose.h"
Chris@82 29 #include <string.h>
Chris@82 30
Chris@82 31 typedef struct {
Chris@82 32 solver super;
Chris@82 33 int preserve_input; /* preserve input even if DESTROY_INPUT was passed */
Chris@82 34 } S;
Chris@82 35
Chris@82 36 typedef struct {
Chris@82 37 plan_mpi_transpose super;
Chris@82 38
Chris@82 39 plan *cld1, *cld2, *cld2rest, *cld3;
Chris@82 40 INT rest_Ioff, rest_Ooff;
Chris@82 41
Chris@82 42 int n_pes, my_pe, *sched;
Chris@82 43 INT *send_block_sizes, *send_block_offsets;
Chris@82 44 INT *recv_block_sizes, *recv_block_offsets;
Chris@82 45 MPI_Comm comm;
Chris@82 46 int preserve_input;
Chris@82 47 } P;
Chris@82 48
Chris@82 49 static void transpose_chunks(int *sched, int n_pes, int my_pe,
Chris@82 50 INT *sbs, INT *sbo, INT *rbs, INT *rbo,
Chris@82 51 MPI_Comm comm,
Chris@82 52 R *I, R *O)
Chris@82 53 {
Chris@82 54 if (sched) {
Chris@82 55 int i;
Chris@82 56 MPI_Status status;
Chris@82 57
Chris@82 58 /* TODO: explore non-synchronous send/recv? */
Chris@82 59
Chris@82 60 if (I == O) {
Chris@82 61 R *buf = (R*) MALLOC(sizeof(R) * sbs[0], BUFFERS);
Chris@82 62
Chris@82 63 for (i = 0; i < n_pes; ++i) {
Chris@82 64 int pe = sched[i];
Chris@82 65 if (my_pe == pe) {
Chris@82 66 if (rbo[pe] != sbo[pe])
Chris@82 67 memmove(O + rbo[pe], O + sbo[pe],
Chris@82 68 sbs[pe] * sizeof(R));
Chris@82 69 }
Chris@82 70 else {
Chris@82 71 memcpy(buf, O + sbo[pe], sbs[pe] * sizeof(R));
Chris@82 72 MPI_Sendrecv(buf, (int) (sbs[pe]), FFTW_MPI_TYPE,
Chris@82 73 pe, (my_pe * n_pes + pe) & 0xffff,
Chris@82 74 O + rbo[pe], (int) (rbs[pe]),
Chris@82 75 FFTW_MPI_TYPE,
Chris@82 76 pe, (pe * n_pes + my_pe) & 0xffff,
Chris@82 77 comm, &status);
Chris@82 78 }
Chris@82 79 }
Chris@82 80
Chris@82 81 X(ifree)(buf);
Chris@82 82 }
Chris@82 83 else { /* I != O */
Chris@82 84 for (i = 0; i < n_pes; ++i) {
Chris@82 85 int pe = sched[i];
Chris@82 86 if (my_pe == pe)
Chris@82 87 memcpy(O + rbo[pe], I + sbo[pe], sbs[pe] * sizeof(R));
Chris@82 88 else
Chris@82 89 MPI_Sendrecv(I + sbo[pe], (int) (sbs[pe]),
Chris@82 90 FFTW_MPI_TYPE,
Chris@82 91 pe, (my_pe * n_pes + pe) & 0xffff,
Chris@82 92 O + rbo[pe], (int) (rbs[pe]),
Chris@82 93 FFTW_MPI_TYPE,
Chris@82 94 pe, (pe * n_pes + my_pe) & 0xffff,
Chris@82 95 comm, &status);
Chris@82 96 }
Chris@82 97 }
Chris@82 98 }
Chris@82 99 }
Chris@82 100
Chris@82 101 static void apply(const plan *ego_, R *I, R *O)
Chris@82 102 {
Chris@82 103 const P *ego = (const P *) ego_;
Chris@82 104 plan_rdft *cld1, *cld2, *cld2rest, *cld3;
Chris@82 105
Chris@82 106 /* transpose locally to get contiguous chunks */
Chris@82 107 cld1 = (plan_rdft *) ego->cld1;
Chris@82 108 if (cld1) {
Chris@82 109 cld1->apply(ego->cld1, I, O);
Chris@82 110
Chris@82 111 if (ego->preserve_input) I = O;
Chris@82 112
Chris@82 113 /* transpose chunks globally */
Chris@82 114 transpose_chunks(ego->sched, ego->n_pes, ego->my_pe,
Chris@82 115 ego->send_block_sizes, ego->send_block_offsets,
Chris@82 116 ego->recv_block_sizes, ego->recv_block_offsets,
Chris@82 117 ego->comm, O, I);
Chris@82 118 }
Chris@82 119 else if (ego->preserve_input) {
Chris@82 120 /* transpose chunks globally */
Chris@82 121 transpose_chunks(ego->sched, ego->n_pes, ego->my_pe,
Chris@82 122 ego->send_block_sizes, ego->send_block_offsets,
Chris@82 123 ego->recv_block_sizes, ego->recv_block_offsets,
Chris@82 124 ego->comm, I, O);
Chris@82 125
Chris@82 126 I = O;
Chris@82 127 }
Chris@82 128 else {
Chris@82 129 /* transpose chunks globally */
Chris@82 130 transpose_chunks(ego->sched, ego->n_pes, ego->my_pe,
Chris@82 131 ego->send_block_sizes, ego->send_block_offsets,
Chris@82 132 ego->recv_block_sizes, ego->recv_block_offsets,
Chris@82 133 ego->comm, I, I);
Chris@82 134 }
Chris@82 135
Chris@82 136 /* transpose locally, again, to get ordinary row-major;
Chris@82 137 this may take two transposes if the block sizes are unequal
Chris@82 138 (3 subplans, two of which operate on disjoint data) */
Chris@82 139 cld2 = (plan_rdft *) ego->cld2;
Chris@82 140 cld2->apply(ego->cld2, I, O);
Chris@82 141 cld2rest = (plan_rdft *) ego->cld2rest;
Chris@82 142 if (cld2rest) {
Chris@82 143 cld2rest->apply(ego->cld2rest,
Chris@82 144 I + ego->rest_Ioff, O + ego->rest_Ooff);
Chris@82 145 cld3 = (plan_rdft *) ego->cld3;
Chris@82 146 if (cld3)
Chris@82 147 cld3->apply(ego->cld3, O, O);
Chris@82 148 /* else TRANSPOSED_OUT is true and user wants O transposed */
Chris@82 149 }
Chris@82 150 }
Chris@82 151
Chris@82 152 static int applicable(const S *ego, const problem *p_,
Chris@82 153 const planner *plnr)
Chris@82 154 {
Chris@82 155 const problem_mpi_transpose *p = (const problem_mpi_transpose *) p_;
Chris@82 156 /* Note: this is *not* UGLY for out-of-place, destroy-input plans;
Chris@82 157 the planner often prefers transpose-pairwise to transpose-alltoall,
Chris@82 158 at least with LAM MPI on my machine. */
Chris@82 159 return (1
Chris@82 160 && (!ego->preserve_input || (!NO_DESTROY_INPUTP(plnr)
Chris@82 161 && p->I != p->O))
Chris@82 162 && ONLY_TRANSPOSEDP(p->flags));
Chris@82 163 }
Chris@82 164
Chris@82 165 static void awake(plan *ego_, enum wakefulness wakefulness)
Chris@82 166 {
Chris@82 167 P *ego = (P *) ego_;
Chris@82 168 X(plan_awake)(ego->cld1, wakefulness);
Chris@82 169 X(plan_awake)(ego->cld2, wakefulness);
Chris@82 170 X(plan_awake)(ego->cld2rest, wakefulness);
Chris@82 171 X(plan_awake)(ego->cld3, wakefulness);
Chris@82 172 }
Chris@82 173
Chris@82 174 static void destroy(plan *ego_)
Chris@82 175 {
Chris@82 176 P *ego = (P *) ego_;
Chris@82 177 X(ifree0)(ego->sched);
Chris@82 178 X(ifree0)(ego->send_block_sizes);
Chris@82 179 MPI_Comm_free(&ego->comm);
Chris@82 180 X(plan_destroy_internal)(ego->cld3);
Chris@82 181 X(plan_destroy_internal)(ego->cld2rest);
Chris@82 182 X(plan_destroy_internal)(ego->cld2);
Chris@82 183 X(plan_destroy_internal)(ego->cld1);
Chris@82 184 }
Chris@82 185
Chris@82 186 static void print(const plan *ego_, printer *p)
Chris@82 187 {
Chris@82 188 const P *ego = (const P *) ego_;
Chris@82 189 p->print(p, "(mpi-transpose-pairwise%s%(%p%)%(%p%)%(%p%)%(%p%))",
Chris@82 190 ego->preserve_input==2 ?"/p":"",
Chris@82 191 ego->cld1, ego->cld2, ego->cld2rest, ego->cld3);
Chris@82 192 }
Chris@82 193
Chris@82 194 /* Given a process which_pe and a number of processes npes, fills
Chris@82 195 the array sched[npes] with a sequence of processes to communicate
Chris@82 196 with for a deadlock-free, optimum-overlap all-to-all communication.
Chris@82 197 (All processes must call this routine to get their own schedules.)
Chris@82 198 The schedule can be re-ordered arbitrarily as long as all processes
Chris@82 199 apply the same permutation to their schedules.
Chris@82 200
Chris@82 201 The algorithm here is based upon the one described in:
Chris@82 202 J. A. M. Schreuder, "Constructing timetables for sport
Chris@82 203 competitions," Mathematical Programming Study 13, pp. 58-67 (1980).
Chris@82 204 In a sport competition, you have N teams and want every team to
Chris@82 205 play every other team in as short a time as possible (maximum overlap
Chris@82 206 between games). This timetabling problem is therefore identical
Chris@82 207 to that of an all-to-all communications problem. In our case, there
Chris@82 208 is one wrinkle: as part of the schedule, the process must do
Chris@82 209 some data transfer with itself (local data movement), analogous
Chris@82 210 to a requirement that each team "play itself" in addition to other
Chris@82 211 teams. With this wrinkle, it turns out that an optimal timetable
Chris@82 212 (N parallel games) can be constructed for any N, not just for even
Chris@82 213 N as in the original problem described by Schreuder.
Chris@82 214 */
Chris@82 215 static void fill1_comm_sched(int *sched, int which_pe, int npes)
Chris@82 216 {
Chris@82 217 int pe, i, n, s = 0;
Chris@82 218 A(which_pe >= 0 && which_pe < npes);
Chris@82 219 if (npes % 2 == 0) {
Chris@82 220 n = npes;
Chris@82 221 sched[s++] = which_pe;
Chris@82 222 }
Chris@82 223 else
Chris@82 224 n = npes + 1;
Chris@82 225 for (pe = 0; pe < n - 1; ++pe) {
Chris@82 226 if (npes % 2 == 0) {
Chris@82 227 if (pe == which_pe) sched[s++] = npes - 1;
Chris@82 228 else if (npes - 1 == which_pe) sched[s++] = pe;
Chris@82 229 }
Chris@82 230 else if (pe == which_pe) sched[s++] = pe;
Chris@82 231
Chris@82 232 if (pe != which_pe && which_pe < n - 1) {
Chris@82 233 i = (pe - which_pe + (n - 1)) % (n - 1);
Chris@82 234 if (i < n/2)
Chris@82 235 sched[s++] = (pe + i) % (n - 1);
Chris@82 236
Chris@82 237 i = (which_pe - pe + (n - 1)) % (n - 1);
Chris@82 238 if (i < n/2)
Chris@82 239 sched[s++] = (pe - i + (n - 1)) % (n - 1);
Chris@82 240 }
Chris@82 241 }
Chris@82 242 A(s == npes);
Chris@82 243 }
Chris@82 244
Chris@82 245 /* Sort the communication schedule sched for npes so that the schedule
Chris@82 246 on process sortpe is ascending or descending (!ascending). This is
Chris@82 247 necessary to allow in-place transposes when the problem does not
Chris@82 248 divide equally among the processes. In this case there is one
Chris@82 249 process where the incoming blocks are bigger/smaller than the
Chris@82 250 outgoing blocks and thus have to be received in
Chris@82 251 descending/ascending order, respectively, to avoid overwriting data
Chris@82 252 before it is sent. */
Chris@82 253 static void sort1_comm_sched(int *sched, int npes, int sortpe, int ascending)
Chris@82 254 {
Chris@82 255 int *sortsched, i;
Chris@82 256 sortsched = (int *) MALLOC(npes * sizeof(int) * 2, OTHER);
Chris@82 257 fill1_comm_sched(sortsched, sortpe, npes);
Chris@82 258 if (ascending)
Chris@82 259 for (i = 0; i < npes; ++i)
Chris@82 260 sortsched[npes + sortsched[i]] = sched[i];
Chris@82 261 else
Chris@82 262 for (i = 0; i < npes; ++i)
Chris@82 263 sortsched[2*npes - 1 - sortsched[i]] = sched[i];
Chris@82 264 for (i = 0; i < npes; ++i)
Chris@82 265 sched[i] = sortsched[npes + i];
Chris@82 266 X(ifree)(sortsched);
Chris@82 267 }
Chris@82 268
Chris@82 269 /* make the plans to do the post-MPI transpositions (shared with
Chris@82 270 transpose-alltoall) */
Chris@82 271 int XM(mkplans_posttranspose)(const problem_mpi_transpose *p, planner *plnr,
Chris@82 272 R *I, R *O, int my_pe,
Chris@82 273 plan **cld2, plan **cld2rest, plan **cld3,
Chris@82 274 INT *rest_Ioff, INT *rest_Ooff)
Chris@82 275 {
Chris@82 276 INT vn = p->vn;
Chris@82 277 INT b = p->block;
Chris@82 278 INT bt = XM(block)(p->ny, p->tblock, my_pe);
Chris@82 279 INT nxb = p->nx / b; /* number of equal-sized blocks */
Chris@82 280 INT nxr = p->nx - nxb * b; /* leftover rows after equal blocks */
Chris@82 281
Chris@82 282 *cld2 = *cld2rest = *cld3 = NULL;
Chris@82 283 *rest_Ioff = *rest_Ooff = 0;
Chris@82 284
Chris@82 285 if (!(p->flags & TRANSPOSED_OUT) && (nxr == 0 || I != O)) {
Chris@82 286 INT nx = p->nx * vn;
Chris@82 287 b *= vn;
Chris@82 288 *cld2 = X(mkplan_f_d)(plnr,
Chris@82 289 X(mkproblem_rdft_0_d)(X(mktensor_3d)
Chris@82 290 (nxb, bt * b, b,
Chris@82 291 bt, b, nx,
Chris@82 292 b, 1, 1),
Chris@82 293 I, O),
Chris@82 294 0, 0, NO_SLOW);
Chris@82 295 if (!*cld2) goto nada;
Chris@82 296
Chris@82 297 if (nxr > 0) {
Chris@82 298 *rest_Ioff = nxb * bt * b;
Chris@82 299 *rest_Ooff = nxb * b;
Chris@82 300 b = nxr * vn;
Chris@82 301 *cld2rest = X(mkplan_f_d)(plnr,
Chris@82 302 X(mkproblem_rdft_0_d)(X(mktensor_2d)
Chris@82 303 (bt, b, nx,
Chris@82 304 b, 1, 1),
Chris@82 305 I + *rest_Ioff,
Chris@82 306 O + *rest_Ooff),
Chris@82 307 0, 0, NO_SLOW);
Chris@82 308 if (!*cld2rest) goto nada;
Chris@82 309 }
Chris@82 310 }
Chris@82 311 else {
Chris@82 312 *cld2 = X(mkplan_f_d)(plnr,
Chris@82 313 X(mkproblem_rdft_0_d)(
Chris@82 314 X(mktensor_4d)
Chris@82 315 (nxb, bt * b * vn, bt * b * vn,
Chris@82 316 bt, b * vn, vn,
Chris@82 317 b, vn, bt * vn,
Chris@82 318 vn, 1, 1),
Chris@82 319 I, O),
Chris@82 320 0, 0, NO_SLOW);
Chris@82 321 if (!*cld2) goto nada;
Chris@82 322
Chris@82 323 *rest_Ioff = *rest_Ooff = nxb * bt * b * vn;
Chris@82 324 *cld2rest = X(mkplan_f_d)(plnr,
Chris@82 325 X(mkproblem_rdft_0_d)(
Chris@82 326 X(mktensor_3d)
Chris@82 327 (bt, nxr * vn, vn,
Chris@82 328 nxr, vn, bt * vn,
Chris@82 329 vn, 1, 1),
Chris@82 330 I + *rest_Ioff, O + *rest_Ooff),
Chris@82 331 0, 0, NO_SLOW);
Chris@82 332 if (!*cld2rest) goto nada;
Chris@82 333
Chris@82 334 if (!(p->flags & TRANSPOSED_OUT)) {
Chris@82 335 *cld3 = X(mkplan_f_d)(plnr,
Chris@82 336 X(mkproblem_rdft_0_d)(
Chris@82 337 X(mktensor_3d)
Chris@82 338 (p->nx, bt * vn, vn,
Chris@82 339 bt, vn, p->nx * vn,
Chris@82 340 vn, 1, 1),
Chris@82 341 O, O),
Chris@82 342 0, 0, NO_SLOW);
Chris@82 343 if (!*cld3) goto nada;
Chris@82 344 }
Chris@82 345 }
Chris@82 346
Chris@82 347 return 1;
Chris@82 348
Chris@82 349 nada:
Chris@82 350 X(plan_destroy_internal)(*cld3);
Chris@82 351 X(plan_destroy_internal)(*cld2rest);
Chris@82 352 X(plan_destroy_internal)(*cld2);
Chris@82 353 *cld2 = *cld2rest = *cld3 = NULL;
Chris@82 354 return 0;
Chris@82 355 }
Chris@82 356
Chris@82 357 static plan *mkplan(const solver *ego_, const problem *p_, planner *plnr)
Chris@82 358 {
Chris@82 359 const S *ego = (const S *) ego_;
Chris@82 360 const problem_mpi_transpose *p;
Chris@82 361 P *pln;
Chris@82 362 plan *cld1 = 0, *cld2 = 0, *cld2rest = 0, *cld3 = 0;
Chris@82 363 INT b, bt, vn, rest_Ioff, rest_Ooff;
Chris@82 364 INT *sbs, *sbo, *rbs, *rbo;
Chris@82 365 int pe, my_pe, n_pes, sort_pe = -1, ascending = 1;
Chris@82 366 R *I, *O;
Chris@82 367 static const plan_adt padt = {
Chris@82 368 XM(transpose_solve), awake, print, destroy
Chris@82 369 };
Chris@82 370
Chris@82 371 UNUSED(ego);
Chris@82 372
Chris@82 373 if (!applicable(ego, p_, plnr))
Chris@82 374 return (plan *) 0;
Chris@82 375
Chris@82 376 p = (const problem_mpi_transpose *) p_;
Chris@82 377 vn = p->vn;
Chris@82 378 I = p->I; O = p->O;
Chris@82 379
Chris@82 380 MPI_Comm_rank(p->comm, &my_pe);
Chris@82 381 MPI_Comm_size(p->comm, &n_pes);
Chris@82 382
Chris@82 383 b = XM(block)(p->nx, p->block, my_pe);
Chris@82 384
Chris@82 385 if (!(p->flags & TRANSPOSED_IN)) { /* b x ny x vn -> ny x b x vn */
Chris@82 386 cld1 = X(mkplan_f_d)(plnr,
Chris@82 387 X(mkproblem_rdft_0_d)(X(mktensor_3d)
Chris@82 388 (b, p->ny * vn, vn,
Chris@82 389 p->ny, vn, b * vn,
Chris@82 390 vn, 1, 1),
Chris@82 391 I, O),
Chris@82 392 0, 0, NO_SLOW);
Chris@82 393 if (XM(any_true)(!cld1, p->comm)) goto nada;
Chris@82 394 }
Chris@82 395 if (ego->preserve_input || NO_DESTROY_INPUTP(plnr)) I = O;
Chris@82 396
Chris@82 397 if (XM(any_true)(!XM(mkplans_posttranspose)(p, plnr, I, O, my_pe,
Chris@82 398 &cld2, &cld2rest, &cld3,
Chris@82 399 &rest_Ioff, &rest_Ooff),
Chris@82 400 p->comm)) goto nada;
Chris@82 401
Chris@82 402 pln = MKPLAN_MPI_TRANSPOSE(P, &padt, apply);
Chris@82 403
Chris@82 404 pln->cld1 = cld1;
Chris@82 405 pln->cld2 = cld2;
Chris@82 406 pln->cld2rest = cld2rest;
Chris@82 407 pln->rest_Ioff = rest_Ioff;
Chris@82 408 pln->rest_Ooff = rest_Ooff;
Chris@82 409 pln->cld3 = cld3;
Chris@82 410 pln->preserve_input = ego->preserve_input ? 2 : NO_DESTROY_INPUTP(plnr);
Chris@82 411
Chris@82 412 MPI_Comm_dup(p->comm, &pln->comm);
Chris@82 413
Chris@82 414 n_pes = (int) X(imax)(XM(num_blocks)(p->nx, p->block),
Chris@82 415 XM(num_blocks)(p->ny, p->tblock));
Chris@82 416
Chris@82 417 /* Compute sizes/offsets of blocks to exchange between processors */
Chris@82 418 sbs = (INT *) MALLOC(4 * n_pes * sizeof(INT), PLANS);
Chris@82 419 sbo = sbs + n_pes;
Chris@82 420 rbs = sbo + n_pes;
Chris@82 421 rbo = rbs + n_pes;
Chris@82 422 b = XM(block)(p->nx, p->block, my_pe);
Chris@82 423 bt = XM(block)(p->ny, p->tblock, my_pe);
Chris@82 424 for (pe = 0; pe < n_pes; ++pe) {
Chris@82 425 INT db, dbt; /* destination block sizes */
Chris@82 426 db = XM(block)(p->nx, p->block, pe);
Chris@82 427 dbt = XM(block)(p->ny, p->tblock, pe);
Chris@82 428
Chris@82 429 sbs[pe] = b * dbt * vn;
Chris@82 430 sbo[pe] = pe * (b * p->tblock) * vn;
Chris@82 431 rbs[pe] = db * bt * vn;
Chris@82 432 rbo[pe] = pe * (p->block * bt) * vn;
Chris@82 433
Chris@82 434 if (db * dbt > 0 && db * p->tblock != p->block * dbt) {
Chris@82 435 A(sort_pe == -1); /* only one process should need sorting */
Chris@82 436 sort_pe = pe;
Chris@82 437 ascending = db * p->tblock > p->block * dbt;
Chris@82 438 }
Chris@82 439 }
Chris@82 440 pln->n_pes = n_pes;
Chris@82 441 pln->my_pe = my_pe;
Chris@82 442 pln->send_block_sizes = sbs;
Chris@82 443 pln->send_block_offsets = sbo;
Chris@82 444 pln->recv_block_sizes = rbs;
Chris@82 445 pln->recv_block_offsets = rbo;
Chris@82 446
Chris@82 447 if (my_pe >= n_pes) {
Chris@82 448 pln->sched = 0; /* this process is not doing anything */
Chris@82 449 }
Chris@82 450 else {
Chris@82 451 pln->sched = (int *) MALLOC(n_pes * sizeof(int), PLANS);
Chris@82 452 fill1_comm_sched(pln->sched, my_pe, n_pes);
Chris@82 453 if (sort_pe >= 0)
Chris@82 454 sort1_comm_sched(pln->sched, n_pes, sort_pe, ascending);
Chris@82 455 }
Chris@82 456
Chris@82 457 X(ops_zero)(&pln->super.super.ops);
Chris@82 458 if (cld1) X(ops_add2)(&cld1->ops, &pln->super.super.ops);
Chris@82 459 if (cld2) X(ops_add2)(&cld2->ops, &pln->super.super.ops);
Chris@82 460 if (cld2rest) X(ops_add2)(&cld2rest->ops, &pln->super.super.ops);
Chris@82 461 if (cld3) X(ops_add2)(&cld3->ops, &pln->super.super.ops);
Chris@82 462 /* FIXME: should MPI operations be counted in "other" somehow? */
Chris@82 463
Chris@82 464 return &(pln->super.super);
Chris@82 465
Chris@82 466 nada:
Chris@82 467 X(plan_destroy_internal)(cld3);
Chris@82 468 X(plan_destroy_internal)(cld2rest);
Chris@82 469 X(plan_destroy_internal)(cld2);
Chris@82 470 X(plan_destroy_internal)(cld1);
Chris@82 471 return (plan *) 0;
Chris@82 472 }
Chris@82 473
Chris@82 474 static solver *mksolver(int preserve_input)
Chris@82 475 {
Chris@82 476 static const solver_adt sadt = { PROBLEM_MPI_TRANSPOSE, mkplan, 0 };
Chris@82 477 S *slv = MKSOLVER(S, &sadt);
Chris@82 478 slv->preserve_input = preserve_input;
Chris@82 479 return &(slv->super);
Chris@82 480 }
Chris@82 481
Chris@82 482 void XM(transpose_pairwise_register)(planner *p)
Chris@82 483 {
Chris@82 484 int preserve_input;
Chris@82 485 for (preserve_input = 0; preserve_input <= 1; ++preserve_input)
Chris@82 486 REGISTER_SOLVER(p, mksolver(preserve_input));
Chris@82 487 }