FORM v5.0.0-35-g6318119
mpidbg.h
Go to the documentation of this file.
1#ifndef MPIDBG_H
2#define MPIDBG_H
3
10/* #[ License : */
11/*
12 * Copyright (C) 1984-2026 J.A.M. Vermaseren
13 * When using this file you are requested to refer to the publication
14 * J.A.M.Vermaseren "New features of FORM" math-ph/0010025
15 * This is considered a matter of courtesy as the development was paid
16 * for by FOM the Dutch physics granting agency and we would like to
17 * be able to track its scientific use to convince FOM of its value
18 * for the community.
19 *
20 * This file is part of FORM.
21 *
22 * FORM is free software: you can redistribute it and/or modify it under the
23 * terms of the GNU General Public License as published by the Free Software
24 * Foundation, either version 3 of the License, or (at your option) any later
25 * version.
26 *
27 * FORM is distributed in the hope that it will be useful, but WITHOUT ANY
28 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
29 * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
30 * details.
31 *
32 * You should have received a copy of the GNU General Public License along
33 * with FORM. If not, see <http://www.gnu.org/licenses/>.
34 */
35/* #] License : */
36
37/*
38 #[ Includes :
39*/
40
41#include <stdarg.h>
42#include <stdio.h>
43#include <string.h>
44#include <mpi.h>
45#if defined(MPIDEBUGGING_DELAY_US) && MPIDEBUGGING_DELAY_US > 0
46#include <unistd.h> // for usleep()
47#endif
48
49#define MPIDBG_LINESIZE 1024
50#define MPIDBG_BUFSIZE 128
51
52/*
53 #] Includes :
54 #[ Utilities :
55 #[ MPIDBG_RANK :
56*/
57
58static inline int MPIDBG_Get_rank(void) {
59 return PF.me; /* Assume we are working with ParFORM. */
60}
61#define MPIDBG_RANK MPIDBG_Get_rank()
62
63/*
64 #] MPIDBG_RANK :
65 #[ MPIDBG_Out :
66*/
67
68static inline void MPIDBG_Out(const char *file, int line, const char *func, const char *fmt, ...) {
69 char buf[MPIDBG_LINESIZE];
70 va_list ap;
71 va_start(ap, fmt);
72 snprintf(buf,MPIDBG_LINESIZE, "*** [%d] %10s %4d @ %-16s: ", MPIDBG_RANK, file, line, func);
73 vsnprintf(buf + strlen(buf),MPIDBG_LINESIZE-strlen(buf), fmt, ap);
74 va_end(ap);
75 /* Assume fprintf with a line will work well even in multi-processes. */
76 fprintf(stderr, "%s\n", buf);
77}
78#define MPIDBG_Out(...) MPIDBG_Out(file, line, func, __VA_ARGS__)
79
80/*
81 #] MPIDBG_Out :
82 #[ MPIDBG_insert_delay :
83*/
84
85static inline void MPIDBG_insert_delay(void)
86{
87#if defined(MPIDEBUGGING_DELAY_US) && MPIDEBUGGING_DELAY_US > 0
88 usleep(MPIDEBUGGING_DELAY_US);
89#endif
90}
91
92/*
93 #] MPIDBG_insert_delay :
94 #[ MPIDBG_sprint_requests :
95*/
96
97static inline void MPIDBG_sprint_requests(char *buf, int count, const MPI_Request *requests)
98{
99 /* Assume sprintf never fail and returns >= 0... */
100 buf += sprintf(buf, "(");
101 int i, first = 1;
102 for ( i = 0; i < count; i++ ) {
103 if ( first ) {
104 first = 0;
105 }
106 else {
107 buf += sprintf(buf, ",");
108 }
109 if ( requests[i] != MPI_REQUEST_NULL ) {
110 buf += sprintf(buf, "%d", i);
111 }
112 else {
113 buf += sprintf(buf, "*");
114 }
115 }
116 buf += sprintf(buf, ")");
117}
118
119/*
120 #] MPIDBG_sprint_requests :
121 #[ MPIDBG_sprint_statuses :
122*/
123
124static inline void MPIDBG_sprint_statuses(char *buf, int count, const MPI_Request *old_requests, const MPI_Request *new_requests, const MPI_Status *statuses)
125{
126 /* Assume sprintf never fail and returns >= 0... */
127 buf += sprintf(buf, "(");
128 int i, first = 1;
129 for ( i = 0; i < count; i++ ) {
130 if ( first ) {
131 first = 0;
132 }
133 else {
134 buf += sprintf(buf, ",");
135 }
136 if ( old_requests[i] != MPI_REQUEST_NULL && new_requests[i] == MPI_REQUEST_NULL ) {
137 int ret_size = 0;
138 MPI_Get_count((MPI_Status *)&statuses[i], MPI_BYTE, &ret_size);
139 buf += sprintf(buf, "(source=%d,tag=%d,size=%d)", statuses[i].MPI_SOURCE, statuses[i].MPI_TAG, ret_size);
140 } else {
141 buf += sprintf(buf, "*");
142 }
143 }
144 buf += sprintf(buf, ")");
145}
146
147/*
148 #] MPIDBG_sprint_statuses :
149 #] Utilities :
150 #[ MPI APIs :
151*/
152
153/*
154 * The followings are the MPI APIs with the logging.
155 */
156
157#define MPIDBG_EXTARG const char *file, int line, const char *func
158
159/*
160 #[ MPI_Init :
161*/
162
163static inline int MPIDBG_Init(int* argc, char*** argv, MPIDBG_EXTARG)
164{
165 int ret = MPI_Init(argc, argv);
166 if ( ret == MPI_SUCCESS ) {
167 MPIDBG_Out("MPI_Init: OK");
168 }
169 else {
170 MPIDBG_Out("MPI_Init: Failed");
171 }
172 return ret;
173}
174#define MPI_Init(...) MPIDBG_Init(__VA_ARGS__, __FILE__, __LINE__, __func__)
175
176/*
177 #] MPI_Init :
178 #[ MPI_Finalize :
179*/
180
181static inline int MPIDBG_Finalize(MPIDBG_EXTARG)
182{
183 MPIDBG_Out("MPI_Finalize");
184 int ret = MPI_Finalize();
185 if ( ret == MPI_SUCCESS ) {
186 MPIDBG_Out("MPI_Finalize: OK");
187 }
188 else {
189 MPIDBG_Out("MPI_Finalize: Failed");
190 }
191 return ret;
192}
193#define MPI_Finalize() MPIDBG_Finalize(__FILE__, __LINE__, __func__)
194
195/*
196 #] MPI_Finalize :
197 #[ MPI_Send :
198*/
199
200static inline int MPIDBG_Send(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
201{
202 MPIDBG_Out("MPI_Send: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
203 int ret = MPI_Send(buf, count, datatype, dest, tag, comm);
204 if ( ret == MPI_SUCCESS ) {
205 MPIDBG_Out("MPI_Send: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
206 }
207 else {
208 MPIDBG_Out("MPI_Send: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
209 }
210 return ret;
211}
212#define MPI_Send(...) MPIDBG_Send(__VA_ARGS__, __FILE__, __LINE__, __func__)
213
214/*
215 #] MPI_Send :
216 #[ MPI_Recv :
217*/
218
219static inline int MPIDBG_Recv(void* buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status *status, MPIDBG_EXTARG)
220{
221 MPI_Status st;
222 if ( status == MPI_STATUS_IGNORE ) status = &st;
223 MPIDBG_Out("MPI_Recv: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
224 int ret = MPI_Recv(buf, count, datatype, source, tag, comm, status);
225 if ( ret == MPI_SUCCESS ) {
226 int ret_count = 0;
227 MPI_Get_count(status, datatype, &ret_count);
228 MPIDBG_Out("MPI_Recv: OK src=%d dest=%d tag=%d count=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_count);
229 }
230 else {
231 MPIDBG_Out("MPI_Recv: Failed src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
232 }
233 return ret;
234}
235#define MPI_Recv(...) MPIDBG_Recv(__VA_ARGS__, __FILE__, __LINE__, __func__)
236
237/*
238 #] MPI_Recv :
239 #[ MPI_Bsend :
240*/
241
242static inline int MPIDBG_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
243{
244 MPIDBG_Out("MPI_Bsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
245 int ret = MPI_Bsend(buf, count, datatype, dest, tag, comm);
246 if ( ret == MPI_SUCCESS ) {
247 MPIDBG_Out("MPI_Bsend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
248 }
249 else {
250 MPIDBG_Out("MPI_Bsend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
251 }
252 return ret;
253}
254#define MPI_Bsend(...) MPIDBG_Bsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
255
256/*
257 #] MPI_Bsend :
258 #[ MPI_Ssend :
259*/
260
261static inline int MPIDBG_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
262{
263 MPIDBG_Out("MPI_Ssend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
264 int ret = MPI_Ssend(buf, count, datatype, dest, tag, comm);
265 if ( ret == MPI_SUCCESS ) {
266 MPIDBG_Out("MPI_Ssend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
267 }
268 else {
269 MPIDBG_Out("MPI_Ssend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
270 }
271 return ret;
272}
273#define MPI_Ssend(...) MPIDBG_Ssend(__VA_ARGS__, __FILE__, __LINE__, __func__)
274
275/*
276 #] MPI_Ssend :
277 #[ MPI_Rsend :
278*/
279
280static inline int MPIDBG_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
281{
282 MPIDBG_Out("MPI_Rsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
283 int ret = MPI_Rsend(buf, count, datatype, dest, tag, comm);
284 if ( ret == MPI_SUCCESS ) {
285 MPIDBG_Out("MPI_Rsend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
286 }
287 else {
288 MPIDBG_Out("MPI_Rsend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
289 }
290 return ret;
291}
292#define MPI_Rsend(...) MPIDBG_Rsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
293
294/*
295 #] MPI_Rsend :
296 #[ MPI_Isend :
297*/
298
299static inline int MPIDBG_Isend(const void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
300{
301 MPIDBG_Out("MPI_Isend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
302 int ret = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
303 if ( ret == MPI_SUCCESS ) {
304 MPIDBG_Out("MPI_Isend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
305 }
306 else {
307 MPIDBG_Out("MPI_Isend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
308 }
309 return ret;
310}
311#define MPI_Isend(...) MPIDBG_Isend(__VA_ARGS__, __FILE__, __LINE__, __func__)
312
313/*
314 #] MPI_Isend :
315 #[ MPI_Ibsend :
316*/
317
318static inline int MPIDBG_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
319{
320 MPIDBG_Out("MPI_Ibsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
321 int ret = MPI_Ibsend(buf, count, datatype, dest, tag, comm, request);
322 if ( ret == MPI_SUCCESS ) {
323 MPIDBG_Out("MPI_Ibsend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
324 }
325 else {
326 MPIDBG_Out("MPI_Ibsend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
327 }
328 return ret;
329}
330#define MPI_Ibsend(...) MPIDBG_Ibsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
331
332/*
333 #] MPI_Ibsend :
334 #[ MPI_Issend :
335*/
336
337static inline int MPIDBG_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
338{
339 MPIDBG_Out("MPI_Issend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
340 int ret = MPI_Issend(buf, count, datatype, dest, tag, comm, request);
341 if ( ret == MPI_SUCCESS ) {
342 MPIDBG_Out("MPI_Issend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
343 }
344 else {
345 MPIDBG_Out("MPI_Issend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
346 }
347 return ret;
348}
349#define MPI_Issend(...) MPIDBG_Issend(__VA_ARGS__, __FILE__, __LINE__, __func__)
350
351/*
352 #] MPI_Issend :
353 #[ MPI_Irsend :
354*/
355
356static inline int MPIDBG_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
357{
358 MPIDBG_Out("MPI_Irsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
359 int ret = MPI_Irsend(buf, count, datatype, dest, tag, comm, request);
360 if ( ret == MPI_SUCCESS ) {
361 MPIDBG_Out("MPI_Irsend: OK src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
362 }
363 else {
364 MPIDBG_Out("MPI_Irsend: Failed src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
365 }
366 return ret;
367}
368#define MPI_Irsend(...) MPIDBG_Irsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
369
370/*
371 #] MPI_Irsend :
372 #[ MPI_Irecv :
373*/
374
375static inline int MPIDBG_Irecv(void* buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
376{
377 MPIDBG_Out("MPI_Irecv: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
378 int ret = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
379 if ( ret == MPI_SUCCESS ) {
380 MPIDBG_Out("MPI_Irecv: OK src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
381 }
382 else {
383 MPIDBG_Out("MPI_Irecv: Failed src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
384 }
385 return ret;
386}
387#define MPI_Irecv(...) MPIDBG_Irecv(__VA_ARGS__, __FILE__, __LINE__, __func__)
388
389/*
390 #] MPI_Irecv :
391 #[ MPI_Wait :
392*/
393
394static inline int MPIDBG_Wait(MPI_Request *request, MPI_Status *status, MPIDBG_EXTARG)
395{
396 MPI_Status st;
397 if ( status == MPI_STATUS_IGNORE ) status = &st;
398 char buf1[MPIDBG_BUFSIZE * 1], buf2[MPIDBG_BUFSIZE * 1];
399 MPI_Request old_request = *request;
400 MPIDBG_sprint_requests(buf1, 1, request);
401 MPIDBG_Out("MPI_Wait: rank=%d request=%s", MPIDBG_RANK, buf1);
402 int ret = MPI_Wait(request, status);
403 if ( ret == MPI_SUCCESS ) {
404 MPIDBG_sprint_statuses(buf2, 1, request, &old_request, status);
405 MPIDBG_Out("MPI_Wait: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
406 }
407 else {
408 MPIDBG_Out("MPI_Wait: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
409 }
410 return ret;
411}
412#define MPI_Wait(...) MPIDBG_Wait(__VA_ARGS__, __FILE__, __LINE__, __func__)
413
414/*
415 #] MPI_Wait :
416 #[ MPI_Test :
417*/
418
419static inline int MPIDBG_Test(MPI_Request *request, int *flag, MPI_Status *status, MPIDBG_EXTARG)
420{
421 MPI_Status st;
422 if ( status == MPI_STATUS_IGNORE ) status = &st;
423 char buf1[MPIDBG_BUFSIZE * 1], buf2[MPIDBG_BUFSIZE * 1];
424 MPI_Request old_request = *request;
425 MPIDBG_sprint_requests(buf1, 1, request);
426 MPIDBG_Out("MPI_Test: rank=%d request=%s", MPIDBG_RANK, buf1);
427 int ret = MPI_Test(request, flag, status);
428 if ( ret == MPI_SUCCESS ) {
429 if ( *flag ) {
430 MPIDBG_sprint_statuses(buf2, 1, request, &old_request, status);
431 MPIDBG_Out("MPI_Test: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
432 }
433 else {
434 MPIDBG_Out("MPI_Test: OK rank=%d request=%s flag=false", MPIDBG_RANK, buf1);
435 }
436 }
437 else {
438 MPIDBG_Out("MPI_Test: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
439 }
440 return ret;
441}
442#define MPI_Test(...) MPIDBG_Test(__VA_ARGS__, __FILE__, __LINE__, __func__)
443
444/*
445 #] MPI_Test :
446 #[ MPI_Waitany :
447*/
448
449static inline int MPIDBG_Waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status, MPIDBG_EXTARG)
450{
451 MPI_Status st;
452 if ( status == MPI_STATUS_IGNORE ) status = &st;
453 char buf1[MPIDBG_BUFSIZE * 1], buf2[MPIDBG_BUFSIZE * 1];
454 MPI_Request old_requests[count];
455 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
456 MPIDBG_sprint_requests(buf1, count, array_of_requests);
457 MPIDBG_Out("MPI_Waitany: rank=%d request=%s", MPIDBG_RANK, buf1);
458 int ret = MPI_Waitany(count, array_of_requests, index, status);
459 if ( ret == MPI_SUCCESS ) {
460 MPI_Status statuses[count];
461 statuses[*index] = *status;
462 MPIDBG_sprint_statuses(buf2, count, old_requests, array_of_requests, statuses);
463 MPIDBG_Out("MPI_Waitany: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
464 }
465 else {
466 MPIDBG_Out("MPI_Waitany: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
467 }
468 return ret;
469}
470#define MPI_Waitany(...) MPIDBG_Waitany(__VA_ARGS__, __FILE__, __LINE__, __func__)
471
472/*
473 #] MPI_Waitany :
474 #[ MPI_Testany :
475*/
476
477static inline int MPIDBG_Testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status, MPIDBG_EXTARG)
478{
479 MPI_Status st;
480 if ( status == MPI_STATUS_IGNORE ) status = &st;
481 char buf1[MPIDBG_BUFSIZE * 1], buf2[MPIDBG_BUFSIZE * 1];
482 MPI_Request old_requests[count];
483 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
484 MPIDBG_sprint_requests(buf1, count, array_of_requests);
485 MPIDBG_Out("MPI_Testany: rank=%d request=%s", MPIDBG_RANK, buf1);
486 int ret = MPI_Testany(count, array_of_requests, index, flag, status);
487 if ( ret == MPI_SUCCESS ) {
488 if ( *flag ) {
489 MPI_Status statuses[count];
490 statuses[*index] = *status;
491 MPIDBG_sprint_statuses(buf2, count, old_requests, array_of_requests, statuses);
492 MPIDBG_Out("MPI_Testany: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
493 }
494 else {
495 MPIDBG_Out("MPI_Testany: OK rank=%d request=%s flag=false", MPIDBG_RANK, buf1);
496 }
497 }
498 else {
499 MPIDBG_Out("MPI_Testany: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
500 }
501 return ret;
502}
503#define MPI_Testany(...) MPIDBG_Testany(__VA_ARGS__, __FILE__, __LINE__, __func__)
504
505/*
506 #] MPI_Testany :
507 #[ MPI_Waitall :
508*/
509
510static inline int MPIDBG_Waitall(int count, MPI_Request *array_of_requests, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
511{
512 MPI_Status st[count];
513 if ( array_of_statuses == MPI_STATUSES_IGNORE ) array_of_statuses = st;
514 char buf1[MPIDBG_BUFSIZE * count], buf2[MPIDBG_BUFSIZE * count];
515 MPI_Request old_requests[count];
516 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
517 MPIDBG_sprint_requests(buf1, count, array_of_requests);
518 MPIDBG_Out("MPI_Waitall: rank=%d request=%s", MPIDBG_RANK, buf1);
519 int ret = MPI_Waitall(count, array_of_requests, array_of_statuses);
520 if ( ret == MPI_SUCCESS ) {
521 MPIDBG_sprint_statuses(buf2, count, old_requests, array_of_requests, array_of_statuses);
522 MPIDBG_Out("MPI_Waitall: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
523 }
524 else {
525 MPIDBG_Out("MPI_Waitall: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
526 }
527 return ret;
528}
529#define MPI_Waitall(...) MPIDBG_Waitall(__VA_ARGS__, __FILE__, __LINE__, __func__)
530
531/*
532 #] MPI_Waitall :
533 #[ MPI_Testall :
534*/
535
536static inline int MPIDBG_Testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
537{
538 MPI_Status st[count];
539 if ( array_of_statuses == MPI_STATUSES_IGNORE ) array_of_statuses = st;
540 char buf1[MPIDBG_BUFSIZE * count], buf2[MPIDBG_BUFSIZE * count];
541 MPI_Request old_requests[count];
542 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
543 MPIDBG_sprint_requests(buf1, count, array_of_requests);
544 MPIDBG_Out("MPI_Testall: rank=%d request=%s", MPIDBG_RANK, buf1);
545 int ret = MPI_Testall(count, array_of_requests, flag, array_of_statuses);
546 if ( ret == MPI_SUCCESS ) {
547 if ( *flag ) {
548 MPIDBG_sprint_statuses(buf2, count, old_requests, array_of_requests, array_of_statuses);
549 MPIDBG_Out("MPI_Testall: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
550 }
551 else {
552 MPIDBG_Out("MPI_Testall: OK rank=%d request=%s flag=false", MPIDBG_RANK, buf1);
553 }
554 }
555 else {
556 MPIDBG_Out("MPI_Testall: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
557 }
558 return ret;
559}
560#define MPI_Testall(...) MPIDBG_Testall(__VA_ARGS__, __FILE__, __LINE__, __func__)
561
562/*
563 #] MPI_Testall :
564 #[ MPI_Waitsome :
565*/
566
567static inline int MPIDBG_Waitsome(int incount, MPI_Request *array_of_requests, int *outcount, int *array_of_indices, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
568{
569 MPI_Status st[incount];
570 if ( array_of_statuses == MPI_STATUSES_IGNORE ) array_of_statuses = st;
571 char buf1[MPIDBG_BUFSIZE * incount], buf2[MPIDBG_BUFSIZE * incount];
572 MPI_Request old_requests[incount];
573 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * incount);
574 MPIDBG_sprint_requests(buf1, incount, array_of_requests);
575 MPIDBG_Out("MPI_Waitsome: rank=%d request=%s", MPIDBG_RANK, buf1);
576 int ret = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
577 if ( ret == MPI_SUCCESS ) {
578 MPIDBG_sprint_statuses(buf2, incount, old_requests, array_of_requests, array_of_statuses);
579 MPIDBG_Out("MPI_Waitsome: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
580 }
581 else {
582 MPIDBG_Out("MPI_Waitsome: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
583 }
584 return ret;
585}
586#define MPI_Waitsome(...) MPIDBG_Waitsome(__VA_ARGS__, __FILE__, __LINE__, __func__)
587
588/*
589 #] MPI_Waitsome :
590 #[ MPI_Testsome :
591*/
592
593static inline int MPIDBG_Testsome(int incount, MPI_Request *array_of_requests, int *outcount, int *array_of_indices, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
594{
595 MPI_Status st[incount];
596 if ( array_of_statuses == MPI_STATUSES_IGNORE ) array_of_statuses = st;
597 char buf1[MPIDBG_BUFSIZE * incount], buf2[MPIDBG_BUFSIZE * incount];
598 MPI_Request old_requests[incount];
599 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * incount);
600 MPIDBG_sprint_requests(buf1, incount, array_of_requests);
601 MPIDBG_Out("MPI_Testsome: rank=%d request=%s", MPIDBG_RANK, buf1);
602 int ret = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
603 if ( ret == MPI_SUCCESS ) {
604 MPIDBG_sprint_statuses(buf2, incount, old_requests, array_of_requests, array_of_statuses);
605 MPIDBG_Out("MPI_Testsome: OK rank=%d request=%s result=%s", MPIDBG_RANK, buf1, buf2);
606 }
607 else {
608 MPIDBG_Out("MPI_Testsome: Failed rank=%d request=%s", MPIDBG_RANK, buf1);
609 }
610 return ret;
611}
612#define MPI_Testsome(...) MPIDBG_Testsome(__VA_ARGS__, __FILE__, __LINE__, __func__)
613
614/*
615 #] MPI_Testsome :
616 #[ MPI_Iprobe :
617*/
618
619static inline int MPIDBG_Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status, MPIDBG_EXTARG)
620{
621 MPI_Status st;
622 if ( status == MPI_STATUS_IGNORE ) status = &st;
623 MPIDBG_Out("MPI_Iprobe: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
624 MPIDBG_insert_delay();
625 int ret = MPI_Iprobe(source, tag, comm, flag, status);
626 if ( ret == MPI_SUCCESS ) {
627 if ( *flag ) {
628 int ret_size = 0;
629 MPI_Get_count(status, MPI_BYTE, &ret_size);
630 MPIDBG_Out("MPI_Iprobe: OK src=%d dest=%d tag=%d size=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_size);
631 }
632 else {
633 MPIDBG_Out("MPI_Iprobe: OK src=%d dest=%d tag=%d flag=false", source, MPIDBG_RANK, tag);
634 }
635 }
636 else {
637 MPIDBG_Out("MPI_Iprobe: Failed src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
638 }
639 return ret;
640}
641#define MPI_Iprobe(...) MPIDBG_Iprobe(__VA_ARGS__, __FILE__, __LINE__, __func__)
642
643/*
644 #] MPI_Iprobe :
645 #[ MPI_Probe :
646*/
647
648static inline int MPIDBG_Probe(int source, int tag, MPI_Comm comm, MPI_Status *status, MPIDBG_EXTARG)
649{
650 MPI_Status st;
651 if ( status == MPI_STATUS_IGNORE ) status = &st;
652 MPIDBG_Out("MPI_Probe: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
653 int ret = MPI_Probe(source, tag, comm, status);
654 if ( ret == MPI_SUCCESS ) {
655 int ret_size = 0;
656 MPI_Get_count(status, MPI_BYTE, &ret_size);
657 MPIDBG_Out("MPI_Probe: OK src=%d dest=%d tag=%d size=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_size);
658 }
659 else {
660 MPIDBG_Out("MPI_Probe: Failed src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
661 }
662 return ret;
663}
664#define MPI_Probe(...) MPIDBG_Probe(__VA_ARGS__, __FILE__, __LINE__, __func__)
665
666/*
667 #] MPI_Probe :
668 #[ MPI_Cancel :
669*/
670
671static inline int MPIDBG_Cancel(MPI_Request *request, MPIDBG_EXTARG)
672{
673 MPIDBG_Out("MPI_Cancel", MPIDBG_RANK);
674 int ret = MPI_Cancel(request);
675 if ( ret == MPI_SUCCESS ) {
676 MPIDBG_Out("MPI_Cancel: OK");
677 }
678 else {
679 MPIDBG_Out("MPI_Cancel: Failed");
680 }
681 return ret;
682}
683#define MPI_Cancel(...) MPIDBG_Cancel(__VA_ARGS__, __FILE__, __LINE__, __func__)
684
685/*
686 #] MPI_Cancel :
687 #[ MPI_Test_cancelled :
688*/
689
690static inline int MPIDBG_Test_cancelled(MPI_Status *status, int *flag, MPIDBG_EXTARG)
691{
692 MPIDBG_Out("MPI_Test_cancelled", MPIDBG_RANK);
693 int ret = MPI_Test_cancelled(status, flag);
694 if ( ret == MPI_SUCCESS ) {
695 if ( *flag ) {
696 MPIDBG_Out("MPI_Test_cancelled: OK flag=true");
697 }
698 else {
699 MPIDBG_Out("MPI_Test_cancelled: OK flag=false");
700 }
701 }
702 else {
703 MPIDBG_Out("MPI_Test_cancelled: Failed");
704 }
705 return ret;
706}
707#define MPI_Test_cancelled(...) MPIDBG_Test_cancelled(__VA_ARGS__, __FILE__, __LINE__, __func__)
708
709/*
710 #] MPI_Test_cancelled :
711 #[ MPI_Barrier :
712*/
713
714static inline int MPIDBG_Barrier(MPI_Comm comm, MPIDBG_EXTARG)
715{
716 MPIDBG_Out("MPI_Barrier");
717 int ret = MPI_Barrier(comm);
718 if ( ret == MPI_SUCCESS ) {
719 MPIDBG_Out("MPI_Barrier: OK");
720 }
721 else {
722 MPIDBG_Out("MPI_Barrier: Failed");
723 }
724 return ret;
725}
726#define MPI_Barrier(...) MPIDBG_Barrier(__VA_ARGS__, __FILE__, __LINE__, __func__)
727
728/*
729 #] MPI_Barrier :
730 #[ MPI_Bcast :
731*/
732
733static inline int MPIDBG_Bcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm, MPIDBG_EXTARG)
734{
735 MPIDBG_Out("MPI_Bcast: root=%d count=%d", root, count);
736 int ret = MPI_Bcast(buffer, count, datatype, root, comm);
737 if ( ret == MPI_SUCCESS ) {
738 MPIDBG_Out("MPI_Bcast: OK root=%d count=%d", root, count);
739 }
740 else {
741 MPIDBG_Out("MPI_Bcast: Failed root=%d count=%d", root, count);
742 }
743 return ret;
744}
745#define MPI_Bcast(...) MPIDBG_Bcast(__VA_ARGS__, __FILE__, __LINE__, __func__)
746
747/*
748 #] MPI_Bcast :
749 #[ MPI_Reduce :
750*/
751
752static inline int MPIDBG_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPIDBG_EXTARG)
753{
754 MPIDBG_Out("MPI_Reduce: root=%d count=%d", root, count);
755 int ret = MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
756 if ( ret == MPI_SUCCESS ) {
757 MPIDBG_Out("MPI_Reduce: OK root=%d count=%d", root, count);
758 }
759 else {
760 MPIDBG_Out("MPI_Reduce: Failed root=%d count=%d", root, count);
761 }
762 return ret;
763}
764#define MPI_Reduce(...) MPIDBG_Reduce(__VA_ARGS__, __FILE__, __LINE__, __func__)
765
766/*
767 #] MPI_Reduce :
768 #] MPI APIs :
769*/
770
771#endif /* MPIDBG_H */