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