annotate fft/fftw/fftw-3.3.4/mpi/wisdom-api.c @ 40:223f770b5341 kissfft-double tip

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