annotate src/fftw-3.3.5/mpi/wisdom-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 7867fa7e1b6b
children
rev   line source
cannam@127 1 /*
cannam@127 2 * Copyright (c) 2003, 2007-14 Matteo Frigo
cannam@127 3 * Copyright (c) 2003, 2007-14 Massachusetts Institute of Technology
cannam@127 4 *
cannam@127 5 * This program is free software; you can redistribute it and/or modify
cannam@127 6 * it under the terms of the GNU General Public License as published by
cannam@127 7 * the Free Software Foundation; either version 2 of the License, or
cannam@127 8 * (at your option) any later version.
cannam@127 9 *
cannam@127 10 * This program is distributed in the hope that it will be useful,
cannam@127 11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
cannam@127 12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
cannam@127 13 * GNU General Public License for more details.
cannam@127 14 *
cannam@127 15 * You should have received a copy of the GNU General Public License
cannam@127 16 * along with this program; if not, write to the Free Software
cannam@127 17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
cannam@127 18 *
cannam@127 19 */
cannam@127 20
cannam@127 21 #include "fftw3-mpi.h"
cannam@127 22 #include "ifftw-mpi.h"
cannam@127 23 #include <string.h>
cannam@127 24
cannam@127 25 #if SIZEOF_SIZE_T == SIZEOF_UNSIGNED_INT
cannam@127 26 # define FFTW_MPI_SIZE_T MPI_UNSIGNED
cannam@127 27 #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG
cannam@127 28 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG
cannam@127 29 #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG_LONG
cannam@127 30 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG_LONG
cannam@127 31 #else
cannam@127 32 # error MPI type for size_t is unknown
cannam@127 33 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG
cannam@127 34 #endif
cannam@127 35
cannam@127 36 /* Import wisdom from all processes to process 0, as prelude to
cannam@127 37 exporting a single wisdom file (this is convenient when we are
cannam@127 38 running on identical processors, to avoid the annoyance of having
cannam@127 39 per-process wisdom files). In order to make the time for this
cannam@127 40 operation logarithmic in the number of processors (rather than
cannam@127 41 linear), we employ a tree reduction algorithm. This means that the
cannam@127 42 wisdom is modified on processes other than root, which shouldn't
cannam@127 43 matter in practice. */
cannam@127 44 void XM(gather_wisdom)(MPI_Comm comm_)
cannam@127 45 {
cannam@127 46 MPI_Comm comm, comm2;
cannam@127 47 int my_pe, n_pes;
cannam@127 48 char *wis;
cannam@127 49 size_t wislen;
cannam@127 50 MPI_Status status;
cannam@127 51
cannam@127 52 MPI_Comm_dup(comm_, &comm);
cannam@127 53 MPI_Comm_rank(comm, &my_pe);
cannam@127 54 MPI_Comm_size(comm, &n_pes);
cannam@127 55
cannam@127 56 if (n_pes > 2) { /* recursively split into even/odd processes */
cannam@127 57 MPI_Comm_split(comm, my_pe % 2, my_pe, &comm2);
cannam@127 58 XM(gather_wisdom)(comm2);
cannam@127 59 MPI_Comm_free(&comm2);
cannam@127 60 }
cannam@127 61 if (n_pes > 1 && my_pe < 2) { /* import process 1 -> 0 */
cannam@127 62 if (my_pe == 1) {
cannam@127 63 wis = X(export_wisdom_to_string)();
cannam@127 64 wislen = strlen(wis) + 1;
cannam@127 65 MPI_Send(&wislen, 1, FFTW_MPI_SIZE_T, 0, 111, comm);
cannam@127 66 MPI_Send(wis, wislen, MPI_CHAR, 0, 222, comm);
cannam@127 67 free(wis);
cannam@127 68 }
cannam@127 69 else /* my_pe == 0 */ {
cannam@127 70 MPI_Recv(&wislen, 1, FFTW_MPI_SIZE_T, 1, 111, comm, &status);
cannam@127 71 wis = (char *) MALLOC(wislen * sizeof(char), OTHER);
cannam@127 72 MPI_Recv(wis, wislen, MPI_CHAR, 1, 222, comm, &status);
cannam@127 73 if (!X(import_wisdom_from_string)(wis))
cannam@127 74 MPI_Abort(comm, 1);
cannam@127 75 X(ifree)(wis);
cannam@127 76 }
cannam@127 77 }
cannam@127 78 MPI_Comm_free(&comm);
cannam@127 79 }
cannam@127 80
cannam@127 81 /* broadcast wisdom from process 0 to all other processes; this
cannam@127 82 is useful so that we can import wisdom once and not worry
cannam@127 83 about parallel I/O or process-specific wisdom, although of
cannam@127 84 course it assumes that all the processes have identical
cannam@127 85 performance characteristics (i.e. identical hardware). */
cannam@127 86 void XM(broadcast_wisdom)(MPI_Comm comm_)
cannam@127 87 {
cannam@127 88 MPI_Comm comm;
cannam@127 89 int my_pe;
cannam@127 90 char *wis;
cannam@127 91 size_t wislen;
cannam@127 92
cannam@127 93 MPI_Comm_dup(comm_, &comm);
cannam@127 94 MPI_Comm_rank(comm, &my_pe);
cannam@127 95
cannam@127 96 if (my_pe != 0) {
cannam@127 97 MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm);
cannam@127 98 wis = (char *) MALLOC(wislen * sizeof(char), OTHER);
cannam@127 99 MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm);
cannam@127 100 if (!X(import_wisdom_from_string)(wis))
cannam@127 101 MPI_Abort(comm, 1);
cannam@127 102 X(ifree)(wis);
cannam@127 103 }
cannam@127 104 else /* my_pe == 0 */ {
cannam@127 105 wis = X(export_wisdom_to_string)();
cannam@127 106 wislen = strlen(wis) + 1;
cannam@127 107 MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm);
cannam@127 108 MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm);
cannam@127 109 X(free)(wis);
cannam@127 110 }
cannam@127 111 MPI_Comm_free(&comm);
cannam@127 112 }