-
Notifications
You must be signed in to change notification settings - Fork 9
/
Copy pathcallee-driven-ex-p1.c
360 lines (342 loc) · 10.7 KB
/
callee-driven-ex-p1.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
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
/*
* Copyright (c) 2013-2019 Triad National Security, LLC
* All rights reserved.
*
* This file is part of the libquo project. See the LICENSE file at the
* top-level directory of this distribution.
*/
#include "callee-driven-ex-p1.h"
#include <stdio.h>
#include <stdlib.h>
#include <stdbool.h>
#include <sys/types.h>
#include <unistd.h>
#include "callee-driven-ex-p1.h"
static int
emit_bind_state(const p1_context_t *c,
char *msg_prefix)
{
char *cbindstr = NULL, *bad_func = NULL;
int bound = 0;
if (QUO_SUCCESS != QUO_stringify_cbind(c->quo, &cbindstr)) {
bad_func = "QUO_stringify_cbind";
goto out;
}
if (QUO_SUCCESS != QUO_bound(c->quo, &bound)) {
bad_func = "QUO_bound";
goto out;
}
printf("%s [rank %d] process %d [%s] bound: %s\n",
msg_prefix, c->rank, (int)getpid(),
cbindstr, bound ? "true" : "false");
fflush(stdout);
out:
if (cbindstr) free(cbindstr);
if (bad_func) {
fprintf(stderr, "%s: %s failure :-(\n", __func__, bad_func);
return 1;
}
return 0;
}
/**
* rudimentary "pretty print" routine. not needed in real life...
*/
static void
demo_emit_sync(const p1_context_t *c)
{
MPI_Barrier(c->init_comm_dup);
usleep((c->rank) * 1000);
}
/**
* gather system and job info from libquo.
*/
static int
sys_grok(p1_context_t *c)
{
char *bad_func = NULL;
/* this interface is more powerful, but the other n* calls can be more
* convenient. at any rate, this is an example of the
* QUO_nobjs_in_type_by_type interface to get the number of sockets on
* the machine. note: you can also use the QUO_nsockets or
* QUO_nobjs_by_type to get the same info. */
if (QUO_SUCCESS != QUO_nobjs_in_type_by_type(c->quo,
QUO_OBJ_MACHINE,
0,
QUO_OBJ_SOCKET,
&c->nsockets)) {
bad_func = "QUO_nobjs_in_type_by_type";
goto out;
}
if (QUO_SUCCESS != QUO_ncores(c->quo, &c->ncores)) {
bad_func = "QUO_ncores";
goto out;
}
if (QUO_SUCCESS != QUO_npus(c->quo, &c->npus)) {
bad_func = "QUO_npus";
goto out;
}
if (QUO_SUCCESS != QUO_bound(c->quo, &c->bound)) {
bad_func = "QUO_bound";
goto out;
}
if (QUO_SUCCESS != QUO_nnodes(c->quo, &c->nnodes)) {
bad_func = "QUO_nnodes";
goto out;
}
if (QUO_SUCCESS != QUO_nqids(c->quo, &c->nnoderanks)) {
bad_func = "QUO_nqids";
goto out;
}
if (QUO_SUCCESS != QUO_id(c->quo, &c->noderank)) {
bad_func = "QUO_id";
goto out;
}
/* for NUMA nodes */
if (QUO_SUCCESS != QUO_nnumanodes(c->quo, &c->nnumanodes)) {
bad_func = "QUO_nnumanodes";
goto out;
}
out:
if (bad_func) {
fprintf(stderr, "%s: %s failure :-(\n", __func__, bad_func);
return 1;
}
return 0;
}
static int
emit_node_basics(const p1_context_t *c)
{
demo_emit_sync(c);
/* one proc per node will emit this info */
if (0 == c->noderank) {
printf("### [rank %d] quo version: %d.%d ###\n",
c->rank, c->qv, c->qsv);
printf("### [rank %d] nnodes: %d\n", c->rank, c->nnodes);
printf("### [rank %d] nnoderanks: %d\n", c->rank, c->nnoderanks);
printf("### [rank %d] nnumanodes: %d\n", c->rank, c->nnumanodes);
printf("### [rank %d] nsockets: %d\n", c->rank, c->nsockets);
printf("### [rank %d] ncores: %d\n", c->rank, c->ncores);
printf("### [rank %d] npus: %d\n", c->rank, c->npus);
fflush(stdout);
}
demo_emit_sync(c);
return 0;
}
/**
* this is where we set our policy regarding who will actually do work. the
* others will sit in a quo barrier an wait for the workers to finish.
*
* this particular example distributes the workers among all the sockets on the
* system, but you can imagine doing the same for NUMA nodes, for example. if
* there are no NUMA nodes on the system, then fall back to something else.
*/
static int
get_worker_pes(p1_context_t *c,
int *nworkers,
int **workers)
{
int res_assigned = 0, tot_workers = 0;
int rc = QUO_ERR;
/* array that holds whether or not a particular rank is going to do work */
int *work_contribs = NULL;
int *worker_ranks = NULL;
/* let quo distribute workers over the sockets. if res_assigned is 1 after
* this call, then i have been chosen. */
if (QUO_SUCCESS != QUO_auto_distrib(c->quo, QUO_OBJ_SOCKET,
2, &res_assigned)) {
return 1;
}
/* array that holds whether or not a particular rank is going to do work */
work_contribs = calloc(c->nranks, sizeof(*work_contribs));
if (!work_contribs) {
rc = QUO_ERR_OOR;
goto out;
}
if (MPI_SUCCESS != (rc = MPI_Allgather(&res_assigned, 1, MPI_INT,
work_contribs, 1, MPI_INT,
c->init_comm_dup))) {
rc = QUO_ERR_MPI;
goto out;
}
/* now iterate over the array and count the total number of workers */
for (int i = 0; i < c->nranks; ++i) {
if (1 == work_contribs[i]) ++tot_workers;
}
worker_ranks = calloc(tot_workers, sizeof(*worker_ranks));
if (!worker_ranks) {
rc = QUO_ERR_OOR;
goto out;
}
/* populate the array with the worker comm world ranks */
for (int i = 0, j = 0; i < c->nranks; ++i) {
if (1 == work_contribs[i]) {
worker_ranks[j++] = i;
}
}
*nworkers = tot_workers;
*workers = worker_ranks;
out:
if (work_contribs) free(work_contribs);
if (QUO_SUCCESS != rc) {
if (worker_ranks) free(worker_ranks);
}
return (QUO_SUCCESS == rc) ? 0 : 1;
}
static int
push_bind(const p1_context_t *c)
{
if (QUO_SUCCESS != QUO_bind_push(c->quo, QUO_BIND_PUSH_OBJ,
QUO_OBJ_SOCKET, -1)) {
return 1;
}
return 0;
}
/* revert our binding policy so p0 can go about its business with its own
* binding policy... */
static int
pop_bind(const p1_context_t *c)
{
if (QUO_SUCCESS != QUO_bind_pop(c->quo)) return 1;
return 0;
}
static int
quo_init(p1_context_t *ctx)
{
/* can be called at any point -- even before construct. */
if (QUO_SUCCESS != QUO_version(&(ctx->qv), &(ctx->qsv))) return 1;
/* relatively expensive call. you only really want to do this once at the
* beginning of time and pass the context all over the place within your
* code.
*/
if (QUO_SUCCESS != QUO_create(&ctx->quo, ctx->init_comm_dup)) return 1;
//
return 0;
}
int
gen_libcomm(p1_context_t *c,
int np1s /* number of participants |p1who| */,
int *p1who /* the participating ranks (initializing comm) */)
{
int rc = QUO_SUCCESS;
if (0 == c->noderank) {
printf("### [rank %d] %d p1pes slated for work\n", c->rank, np1s);
printf("### [rank %d] and they are: ", c->rank);
if (0 == np1s) printf("\n");
fflush(stdout);
for (int i = 0; i < np1s; ++i) {
printf("%d ", p1who[i]);
fflush(stdout);
if (i + 1 == np1s) {
printf("\n");
fflush(stdout);
}
}
}
/* ////////////////////////////////////////////////////////////////////// */
/* now create our own communicator based on the rank ids passed here */
/* ////////////////////////////////////////////////////////////////////// */
MPI_Group init_comm_grp;
MPI_Group p1_group;
if (MPI_SUCCESS != MPI_Comm_group(c->init_comm_dup, &init_comm_grp)) {
rc = QUO_ERR_MPI;
goto out;
}
if (MPI_SUCCESS != MPI_Group_incl(init_comm_grp, np1s,
p1who, &p1_group)) {
rc = QUO_ERR_MPI;
goto out;
}
if (MPI_SUCCESS != MPI_Comm_create(c->init_comm_dup,
p1_group,
&(c->quo_comm))) {
rc = QUO_ERR_MPI;
goto out;
}
/* am i in the new communicator? */
c->in_quo_comm = (MPI_COMM_NULL == c->quo_comm) ? false : true;
if (c->in_quo_comm) {
if (MPI_SUCCESS != MPI_Comm_size(c->quo_comm, &c->qc_size)) {
rc = QUO_ERR_MPI;
goto out;
}
if (MPI_SUCCESS != MPI_Comm_rank(c->quo_comm, &c->qc_rank)) {
rc = QUO_ERR_MPI;
goto out;
}
}
out:
if (MPI_SUCCESS != MPI_Group_free(&init_comm_grp)) return 1;
return (QUO_SUCCESS == rc) ? 0 : 1;
}
int
p1_init(p1_context_t **p1_ctx,
MPI_Comm comm)
{
if (!p1_ctx) return 1;
//
int inited = 0;
if (MPI_SUCCESS != MPI_Initialized(&inited)) return 1;
/* QUO requires that MPI be initialized before its use. */
if (!inited) return 1;
//
p1_context_t *newc = NULL;
if (NULL == (newc = calloc(1, sizeof(*newc)))) return 1;
// dup initializing comm */
if (MPI_SUCCESS != MPI_Comm_dup(comm, &newc->init_comm_dup)) return 1;
/* gather some basic info about initializing communicator */
if (MPI_SUCCESS != MPI_Comm_size(newc->init_comm_dup, &newc->nranks)) {
return 1;
}
if (MPI_SUCCESS != MPI_Comm_rank(newc->init_comm_dup, &newc->rank)) {
return 1;
}
//
if (quo_init(newc)) return 1;
//
if (sys_grok(newc)) return 1;
//
if (emit_node_basics(newc)) return 1;
//
demo_emit_sync(newc);
if (emit_bind_state(newc, "###")) return 1;
demo_emit_sync(newc);
//
int n_workers = 0;
int *worker_comm_ids = NULL;
if (get_worker_pes(newc, &n_workers, &worker_comm_ids)) goto err;
if (gen_libcomm(newc, n_workers, worker_comm_ids)) goto err;
//
*p1_ctx = newc;
if (worker_comm_ids) free(worker_comm_ids);
return 0;
err:
if (worker_comm_ids) free(worker_comm_ids);
return 1;
}
int
p1_fini(p1_context_t *c)
{
if (!c) return 0;
MPI_Comm_free(&c->init_comm_dup);
if (c->in_quo_comm) MPI_Comm_free(&c->quo_comm);
QUO_free(c->quo);
free(c);
return 0;
}
int
p1_entry_point(p1_context_t *c)
{
/* actually do threaded work? */
if (c->in_quo_comm) {
/* prep runtime environment */
push_bind(c);
/* DO WORK HERE */
if (emit_bind_state(c, "-->")) return 1;
/* revert to previous runtime environment */
pop_bind(c);
}
if (QUO_SUCCESS != QUO_barrier(c->quo)) return 1;
/* via mpi, share results with processes that were not in quo_comm, or any
* other processes that need those data. sharing is caring. */
return 0;
}