-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathcomm.c
345 lines (309 loc) · 9.09 KB
/
comm.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
/*
This file is part of SystemConfidence.
Copyright (C) 2012, UT-Battelle, LLC.
This product includes software produced by UT-Battelle, LLC under Contract No.
DE-AC05-00OR22725 with the Department of Energy.
This program is free software; you can redistribute it and/or modify
it under the terms of the New BSD 3-clause software license (LICENSE).
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
LICENSE for more details.
For more information please contact the SystemConfidence developers at:
systemconfidence-info@googlegroups.com
*/
#include <unistd.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <stdio.h>
#include <errno.h>
#include <ctype.h>
#include "config.h"
#include "orbtimer.h"
#include "comm.h"
#ifdef SHMEM
#include <mpp/shmem.h>
# ifndef _SHMEM_COLLECT_SYNC_SIZE
// FIXME: hacks to cope with Cray....
# define _SHMEM_COLLECT_SYNC_SIZE _SHMEM_REDUCE_SYNC_SIZE
# define _my_pe my_pe
# define _num_pes num_pes
# endif
static long cSync[_SHMEM_COLLECT_SYNC_SIZE];
static long rSync[_SHMEM_REDUCE_SYNC_SIZE];
#else
#include <mpi.h>
#endif
/**
* \brief routines in this file abstract the communication layer.
*
* so far, only MPI is implemented
*
*/
/**
* \brief Creates a buffer struct to hold the communication info
* \param nbytes Tells the size of the data array
* \return p Buffer variable
*/
buffer_p comm_newbuffer(size_t nbytes) {
buffer_p p;
p = (buffer_p) malloc(sizeof(buffer_t));
assert(p != NULL);
p->len = nbytes;
#ifdef SHMEM
p->data = (void *)shmalloc(nbytes);
#else /* MPI */
p->data = (void *)malloc(nbytes);
#endif
assert(p->data != NULL);
return p;
}
/**
* \brief Frees the buffer struct
* \param buf Buffer to be free'd
* \sa comm_newbuffer
*/
void comm_freebuffer(buffer_p buf) {
#ifdef SHMEM
shfree(buf->data);
#else /* MPI */
free(buf->data);
#endif
free(buf);
return;
}
/**
* \brief Allocates space for the dist array
* \param num_bins The size of the dist array
*/
uint64_t *comm_alloc_dist(size_t num_bins) {
uint64_t *dist_array;
#ifdef SHMEM
int i;
dist_array = (uint64_t *) shmalloc(num_bins * sizeof(uint64_t));
assert(dist_array != NULL);
for (i = 0; i < num_bins; i++) {
dist_array[i] = 0x0;
}
#else /* MPI */
dist_array = (uint64_t *) calloc(num_bins, sizeof(uint64_t));
#endif
return dist_array;
}
/**
* \brief Frees the dist array
* \param dist_array Memory location to be free'd
* \sa comm_alloc_dist
*/
void comm_free_dist(uint64_t *dist_array) {
#ifdef SHMEM
shfree(dist_array);
#else /* MPI */
free(dist_array);
#endif
}
/**
* \brief Collects the measurements into a global location
* \param g The global array of measurements collected over the course of the run
* \param l The local array of measurements
*/
void comm_aggregate(measurement_p g, measurement_p l) {
/* collects local measurments into a global measurement */
int i, ierr;
assert(l->nbins == g->nbins);
assert(l->num_histograms == g->num_histograms);
ierr = 0;
#ifdef SHMEM
int max = (l->nbins/2 + 1) > _SHMEM_REDUCE_MIN_WRKDATA_SIZE ? (l->nbins/2 + 1) : _SHMEM_REDUCE_MIN_WRKDATA_SIZE;
uint64_t *pWrk = (uint64_t *)shmalloc(max * sizeof(uint64_t));
assert(pWrk != NULL);
for (i = 0; i < l->num_histograms; i++) {
shmem_barrier_all();
shmem_longlong_sum_to_all((long long *)g->hist[i].dist, (long long *)l->hist[i].dist,
l->nbins, 0, 0, num_ranks, (long long *)pWrk, rSync);
}
shmem_barrier_all();
shfree(pWrk);
#else /* MPI case */
for (i = 0; i < l->num_histograms; i++) {
ierr += MPI_Allreduce(l->hist[i].dist, g->hist[i].dist, l->nbins,
MPI_INTEGER8, MPI_SUM, MPI_COMM_WORLD);
}
assert(ierr == 0);
#endif
return;
}
/**
* \brief Prints out the mapping of the test results
* \param tst Holds the information used for printing
*/
void comm_showmapping(test_p tst) {
int i;
FILE *Fmapping;
char fname[512];
ROOTONLY {
if (tst->rank_mapping == 1) {
snprintf(fname, 512, "%s/%s.RankToNodeMapping.%dof%d",
tst->case_name, "global", my_rank, num_ranks);
Fmapping = fopen(fname, "w");
fprintf(Fmapping, "#Task Rank, node_ids, (from RootTask)\n");
for (i = 0; i < num_ranks; i++)
fprintf(Fmapping,"%d, %ld\n", i, node_id[i]);
fclose(Fmapping);
}
}
return;
}
/**
* \brief Gives the value of the calling node
* \return id Is the ID of the node
*/
uint64_t comm_getnodeid() {
/*
* Assuming the system's nodes are "numbered" after some fashion, this
* function returns a unique integer for each hardware 'node'. The
* nodeid's of two MPI ranks are compared and the communication is
* assumed to be local (ie. at a latency potentially less than that of
* the network) if the nodeid's are the same. If something other than
* digits are used to differentiate the nodes this function should be
* replaced with a function which generates a useful nodeid for the
* target machine. If none of the NODEID_<TYPE> macros are set, we
* punt, and return nodeid's suggesting every rank on a separate node.
*/
uint64_t id;
#if defined(NODEID_GETHOSTNAME) || defined(NODEID_MPI) || defined(NODEID_SLURM)
char *pn, *pp, *limit;
int ierr = 0;
# if defined(SHMEM)
gethostname(namebuff, NAMEBUFFSIZE);
nodename = namebuff;
# elif defined(NODEID_MPI)
int len = NAMEBUFFSIZE;
ierr = MPI_Get_processor_name(namebuff, &len);
nodename = namebuff;
# elif defined(NODEID_GETHOSTNAME)
gethostname(namebuff, NAMEBUFFSIZE);
nodename = namebuff;
# elif defined(NODEID_SLURM)
nodename = getenv("SLURM_NODEID");
# endif
limit = nodename + strlen(nodename);
pp = nodename - 1;
pn = nid;
/* copy all (and only) the digits to a new buffer */
while (++pp < nodename + strlen(nodename))
if (isdigit(*pp))
*pn++ = *pp;
*pn = (char)0;
errno = 0;
id = strtoull(nid, NULL, 10);
ierr = errno;
if (ierr == ERANGE || ierr == EINVAL) {
perror("strtol");
exit(ierr);
}
#else /* NODEID_<TYPE> not defined. Punt. One MPI rank per node. */
int tmp;
# ifdef SHMEM
tmp = _my_pe();
# else
int ierr = 0;
ierr = MPI_Comm_rank(MPI_COMM_WORLD, &tmp);
assert(ierr == 0);
# endif
id = tmp;
#endif /* NODEID_<TYPE> */
return id;
}
/**
* \brief Returns the next power of 2 greater than n
* \param n 2^c > n, 2^(c-1) <= n
* \return c 2^c > n, 2^(c-1) <= n
*/
int comm_ceil2(int n) {
int c = 1;
while (c < n)
c *= 2;
return c;
}
/**
* \brief Initializes the communication arrays - SHMEM
* \param tst Will tell the test how many stages it should run
* \param argc Typical C variable: tells the number of arguments in the command line
* \param argv Typical C variable: holds the command line arguments
*/
void comm_SHMEM_initialize(test_p tst, int *argc, char **argv[]) {
#ifdef SHMEM
uint64_t mynodeid;
int i;
start_pes(0);
my_rank = _my_pe();
num_ranks = _num_pes();
root_rank = 0;
tst->num_stages = comm_ceil2(num_ranks);
node_id = (uint64_t *)shmalloc(num_ranks * sizeof(uint64_t));
assert(node_id != NULL);
for ( i=0; i<num_ranks; i++) node_id[i] = 0;
/* sync arrays for collectives */
for (i = 0; i < _SHMEM_COLLECT_SYNC_SIZE; i++) cSync[i] = _SHMEM_SYNC_VALUE;
for (i = 0; i < _SHMEM_REDUCE_SYNC_SIZE; i++) rSync[i] = _SHMEM_SYNC_VALUE;
shmem_barrier_all();
node_id[my_rank] = mynodeid = comm_getnodeid();
// FIXME:
// would prefer collect64, but it doesn't exist on some platforms. is shmem_longlong_sum_to_all() better?
// shmem_collect64((void *)node_id, (void *)&mynodeid, 1, 0, 0, num_ranks, cSync);
int max = (1+num_ranks/2)>_SHMEM_REDUCE_MIN_WRKDATA_SIZE ? (1+num_ranks/2) : _SHMEM_REDUCE_MIN_WRKDATA_SIZE;
uint64_t *pWrk = (uint64_t *)shmalloc(max * sizeof(uint64_t));
assert(pWrk != NULL);
shmem_longlong_sum_to_all((long long *)node_id, (long long *)node_id, num_ranks, 0, 0, num_ranks, pWrk, rSync);
shmem_barrier_all();
#endif /* SHMEM */
return;
}
/**
* \brief Frees up the space used by comm - SHMEM
* \sa comm_SHMEM_initialize
*/
void comm_SHMEM_finalize() {
#ifdef SHMEM
shfree(node_id);
#endif
/* SHMEM requires no finalize function */
return;
}
/**
* \brief Initializes the communication arrays for comm - MPI
* \param tst Will tell the test how many stages it should run
* \param argc Typical C variable: tells the number of arguments in the command line
* \param argv Typical C variable: holds the command line arguments
*/
void comm_MPI_initialize(test_p tst, int *argc, char **argv[]) {
#ifndef SHMEM
uint64_t mynodeid;
int ierr = 0;
ierr += MPI_Init(argc, argv);
ierr += MPI_Comm_rank(MPI_COMM_WORLD, &my_rank);
ierr += MPI_Comm_size(MPI_COMM_WORLD, &num_ranks);
root_rank = 0;
tst->num_stages = comm_ceil2(num_ranks);
node_id = (uint64_t *)malloc(num_ranks * sizeof(uint64_t));
assert(node_id != NULL);
node_id[my_rank] = mynodeid = comm_getnodeid();
ierr += MPI_Allgather(&mynodeid, sizeof(uint64_t), MPI_BYTE,
node_id, sizeof(uint64_t), MPI_BYTE, MPI_COMM_WORLD);
assert(ierr == 0);
#endif
return;
}
/**
* \brief Frees up the spave used by comm - MPI
* \sa comm_MPI_initialize
*/
void comm_MPI_finalize() {
#ifndef SHMEM
free(node_id);
MPI_Finalize();
#endif
return;
}