Chris@10: /* Chris@10: * Copyright (c) 2003, 2007-11 Matteo Frigo Chris@10: * Copyright (c) 2003, 2007-11 Massachusetts Institute of Technology Chris@10: * Chris@10: * This program is free software; you can redistribute it and/or modify Chris@10: * it under the terms of the GNU General Public License as published by Chris@10: * the Free Software Foundation; either version 2 of the License, or Chris@10: * (at your option) any later version. Chris@10: * Chris@10: * This program is distributed in the hope that it will be useful, Chris@10: * but WITHOUT ANY WARRANTY; without even the implied warranty of Chris@10: * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the Chris@10: * GNU General Public License for more details. Chris@10: * Chris@10: * You should have received a copy of the GNU General Public License Chris@10: * along with this program; if not, write to the Free Software Chris@10: * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA Chris@10: * Chris@10: */ Chris@10: Chris@10: #include "fftw3-mpi.h" Chris@10: #include "ifftw-mpi.h" Chris@10: #include Chris@10: Chris@10: #if SIZEOF_SIZE_T == SIZEOF_UNSIGNED_INT Chris@10: # define FFTW_MPI_SIZE_T MPI_UNSIGNED Chris@10: #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG Chris@10: # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG Chris@10: #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG_LONG Chris@10: # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG_LONG Chris@10: #else Chris@10: # error MPI type for size_t is unknown Chris@10: # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG Chris@10: #endif Chris@10: Chris@10: /* Import wisdom from all processes to process 0, as prelude to Chris@10: exporting a single wisdom file (this is convenient when we are Chris@10: running on identical processors, to avoid the annoyance of having Chris@10: per-process wisdom files). In order to make the time for this Chris@10: operation logarithmic in the number of processors (rather than Chris@10: linear), we employ a tree reduction algorithm. This means that the Chris@10: wisdom is modified on processes other than root, which shouldn't Chris@10: matter in practice. */ Chris@10: void XM(gather_wisdom)(MPI_Comm comm_) Chris@10: { Chris@10: MPI_Comm comm, comm2; Chris@10: int my_pe, n_pes; Chris@10: char *wis; Chris@10: size_t wislen; Chris@10: MPI_Status status; Chris@10: Chris@10: MPI_Comm_dup(comm_, &comm); Chris@10: MPI_Comm_rank(comm, &my_pe); Chris@10: MPI_Comm_size(comm, &n_pes); Chris@10: Chris@10: if (n_pes > 2) { /* recursively split into even/odd processes */ Chris@10: MPI_Comm_split(comm, my_pe % 2, my_pe, &comm2); Chris@10: XM(gather_wisdom)(comm2); Chris@10: MPI_Comm_free(&comm2); Chris@10: } Chris@10: if (n_pes > 1 && my_pe < 2) { /* import process 1 -> 0 */ Chris@10: if (my_pe == 1) { Chris@10: wis = X(export_wisdom_to_string)(); Chris@10: wislen = strlen(wis) + 1; Chris@10: MPI_Send(&wislen, 1, FFTW_MPI_SIZE_T, 0, 111, comm); Chris@10: MPI_Send(wis, wislen, MPI_CHAR, 0, 222, comm); Chris@10: free(wis); Chris@10: } Chris@10: else /* my_pe == 0 */ { Chris@10: MPI_Recv(&wislen, 1, FFTW_MPI_SIZE_T, 1, 111, comm, &status); Chris@10: wis = (char *) MALLOC(wislen * sizeof(char), OTHER); Chris@10: MPI_Recv(wis, wislen, MPI_CHAR, 1, 222, comm, &status); Chris@10: if (!X(import_wisdom_from_string)(wis)) Chris@10: MPI_Abort(comm, 1); Chris@10: X(ifree)(wis); Chris@10: } Chris@10: } Chris@10: MPI_Comm_free(&comm); Chris@10: } Chris@10: Chris@10: /* broadcast wisdom from process 0 to all other processes; this Chris@10: is useful so that we can import wisdom once and not worry Chris@10: about parallel I/O or process-specific wisdom, although of Chris@10: course it assumes that all the processes have identical Chris@10: performance characteristics (i.e. identical hardware). */ Chris@10: void XM(broadcast_wisdom)(MPI_Comm comm_) Chris@10: { Chris@10: MPI_Comm comm; Chris@10: int my_pe; Chris@10: char *wis; Chris@10: size_t wislen; Chris@10: Chris@10: MPI_Comm_dup(comm_, &comm); Chris@10: MPI_Comm_rank(comm, &my_pe); Chris@10: Chris@10: if (my_pe != 0) { Chris@10: MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm); Chris@10: wis = (char *) MALLOC(wislen * sizeof(char), OTHER); Chris@10: MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm); Chris@10: if (!X(import_wisdom_from_string)(wis)) Chris@10: MPI_Abort(comm, 1); Chris@10: X(ifree)(wis); Chris@10: } Chris@10: else /* my_pe == 0 */ { Chris@10: wis = X(export_wisdom_to_string)(); Chris@10: wislen = strlen(wis) + 1; Chris@10: MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm); Chris@10: MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm); Chris@10: X(free)(wis); Chris@10: } Chris@10: MPI_Comm_free(&comm); Chris@10: }