annotate src/fftw-3.3.3/mpi/wisdom-api.c @ 95:89f5e221ed7b

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