annotate src/fftw-3.3.8/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 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 "fftw3-mpi.h"
cannam@167 22 #include "ifftw-mpi.h"
cannam@167 23 #include <string.h>
cannam@167 24
cannam@167 25 #if SIZEOF_SIZE_T == SIZEOF_UNSIGNED_INT
cannam@167 26 # define FFTW_MPI_SIZE_T MPI_UNSIGNED
cannam@167 27 #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG
cannam@167 28 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG
cannam@167 29 #elif SIZEOF_SIZE_T == SIZEOF_UNSIGNED_LONG_LONG
cannam@167 30 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG_LONG
cannam@167 31 #else
cannam@167 32 # error MPI type for size_t is unknown
cannam@167 33 # define FFTW_MPI_SIZE_T MPI_UNSIGNED_LONG
cannam@167 34 #endif
cannam@167 35
cannam@167 36 /* Import wisdom from all processes to process 0, as prelude to
cannam@167 37 exporting a single wisdom file (this is convenient when we are
cannam@167 38 running on identical processors, to avoid the annoyance of having
cannam@167 39 per-process wisdom files). In order to make the time for this
cannam@167 40 operation logarithmic in the number of processors (rather than
cannam@167 41 linear), we employ a tree reduction algorithm. This means that the
cannam@167 42 wisdom is modified on processes other than root, which shouldn't
cannam@167 43 matter in practice. */
cannam@167 44 void XM(gather_wisdom)(MPI_Comm comm_)
cannam@167 45 {
cannam@167 46 MPI_Comm comm, comm2;
cannam@167 47 int my_pe, n_pes;
cannam@167 48 char *wis;
cannam@167 49 size_t wislen;
cannam@167 50 MPI_Status status;
cannam@167 51
cannam@167 52 MPI_Comm_dup(comm_, &comm);
cannam@167 53 MPI_Comm_rank(comm, &my_pe);
cannam@167 54 MPI_Comm_size(comm, &n_pes);
cannam@167 55
cannam@167 56 if (n_pes > 2) { /* recursively split into even/odd processes */
cannam@167 57 MPI_Comm_split(comm, my_pe % 2, my_pe, &comm2);
cannam@167 58 XM(gather_wisdom)(comm2);
cannam@167 59 MPI_Comm_free(&comm2);
cannam@167 60 }
cannam@167 61 if (n_pes > 1 && my_pe < 2) { /* import process 1 -> 0 */
cannam@167 62 if (my_pe == 1) {
cannam@167 63 wis = X(export_wisdom_to_string)();
cannam@167 64 wislen = strlen(wis) + 1;
cannam@167 65 MPI_Send(&wislen, 1, FFTW_MPI_SIZE_T, 0, 111, comm);
cannam@167 66 MPI_Send(wis, wislen, MPI_CHAR, 0, 222, comm);
cannam@167 67 free(wis);
cannam@167 68 }
cannam@167 69 else /* my_pe == 0 */ {
cannam@167 70 MPI_Recv(&wislen, 1, FFTW_MPI_SIZE_T, 1, 111, comm, &status);
cannam@167 71 wis = (char *) MALLOC(wislen * sizeof(char), OTHER);
cannam@167 72 MPI_Recv(wis, wislen, MPI_CHAR, 1, 222, comm, &status);
cannam@167 73 if (!X(import_wisdom_from_string)(wis))
cannam@167 74 MPI_Abort(comm, 1);
cannam@167 75 X(ifree)(wis);
cannam@167 76 }
cannam@167 77 }
cannam@167 78 MPI_Comm_free(&comm);
cannam@167 79 }
cannam@167 80
cannam@167 81 /* broadcast wisdom from process 0 to all other processes; this
cannam@167 82 is useful so that we can import wisdom once and not worry
cannam@167 83 about parallel I/O or process-specific wisdom, although of
cannam@167 84 course it assumes that all the processes have identical
cannam@167 85 performance characteristics (i.e. identical hardware). */
cannam@167 86 void XM(broadcast_wisdom)(MPI_Comm comm_)
cannam@167 87 {
cannam@167 88 MPI_Comm comm;
cannam@167 89 int my_pe;
cannam@167 90 char *wis;
cannam@167 91 size_t wislen;
cannam@167 92
cannam@167 93 MPI_Comm_dup(comm_, &comm);
cannam@167 94 MPI_Comm_rank(comm, &my_pe);
cannam@167 95
cannam@167 96 if (my_pe != 0) {
cannam@167 97 MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm);
cannam@167 98 wis = (char *) MALLOC(wislen * sizeof(char), OTHER);
cannam@167 99 MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm);
cannam@167 100 if (!X(import_wisdom_from_string)(wis))
cannam@167 101 MPI_Abort(comm, 1);
cannam@167 102 X(ifree)(wis);
cannam@167 103 }
cannam@167 104 else /* my_pe == 0 */ {
cannam@167 105 wis = X(export_wisdom_to_string)();
cannam@167 106 wislen = strlen(wis) + 1;
cannam@167 107 MPI_Bcast(&wislen, 1, FFTW_MPI_SIZE_T, 0, comm);
cannam@167 108 MPI_Bcast(wis, wislen, MPI_CHAR, 0, comm);
cannam@167 109 X(free)(wis);
cannam@167 110 }
cannam@167 111 MPI_Comm_free(&comm);
cannam@167 112 }