FORM 4.3
mpidbg.h
Go to the documentation of this file.
1#ifndef MPIDBG_H
2#define MPIDBG_H
3
9
10/* #[ License : */
11/*
12 * Copyright (C) 1984-2022 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
46/*
47 #] Includes :
48 #[ Utilities :
49 #[ MPIDBG_RANK :
50*/
51
52static inline int MPIDBG_Get_rank(void) {
53 return PF.me; /* Assume we are working with ParFORM. */
54}
55#define MPIDBG_RANK MPIDBG_Get_rank()
56
57/*
58 #] MPIDBG_RANK :
59 #[ MPIDBG_Out :
60*/
61
62static inline void MPIDBG_Out(const char *file, int line, const char *func, const char *fmt, ...) {
63 char buf[1024]; /* Enough. */
64 va_list ap;
65 va_start(ap, fmt);
66 sprintf(buf, "*** [%d] %10s %4d @ %-16s: ", MPIDBG_RANK, file, line, func);
67 vsprintf(buf + strlen(buf), fmt, ap);
68 va_end(ap);
69 /* Assume fprintf with a line will work well even in multi-processes. */
70 fprintf(stderr, "%s\n", buf);
71}
72#define MPIDBG_Out(...) MPIDBG_Out(file, line, func, __VA_ARGS__)
73
74/*
75 #] MPIDBG_Out :
76 #[ MPIDBG_sprint_requests :
77*/
78
79static inline void MPIDBG_sprint_requests(char *buf, int count, const MPI_Request *requests)
80{
81 /* Assume sprintf never fail and returns >= 0... */
82 buf += sprintf(buf, "(");
83 int i, first = 1;
84 for ( i = 0; i < count; i++ ) {
85 if ( first ) {
86 first = 0;
87 }
88 else {
89 buf += sprintf(buf, ",");
90 }
91 if ( requests[i] != MPI_REQUEST_NULL ) {
92 buf += sprintf(buf, "%d", i);
93 }
94 else {
95 buf += sprintf(buf, "*");
96 }
97 }
98 buf += sprintf(buf, ")");
99}
100
101/*
102 #] MPIDBG_sprint_requests :
103 #[ MPIDBG_sprint_statuses :
104*/
105
106static inline void MPIDBG_sprint_statuses(char *buf, int count, const MPI_Request *old_requests, const MPI_Request *new_requests, const MPI_Status *statuses)
107{
108 /* Assume sprintf never fail and returns >= 0... */
109 buf += sprintf(buf, "(");
110 int i, first = 1;
111 for ( i = 0; i < count; i++ ) {
112 if ( first ) {
113 first = 0;
114 }
115 else {
116 buf += sprintf(buf, ",");
117 }
118 if ( old_requests[i] != MPI_REQUEST_NULL && new_requests[i] == MPI_REQUEST_NULL ) {
119 int ret_size = 0;
120 MPI_Get_count((MPI_Status *)&statuses[i], MPI_BYTE, &ret_size);
121 buf += sprintf(buf, "(source=%d,tag=%d,size=%d)", statuses[i].MPI_SOURCE, statuses[i].MPI_TAG, ret_size);
122 } else {
123 buf += sprintf(buf, "*");
124 }
125 }
126 buf += sprintf(buf, ")");
127}
128
129/*
130 #] MPIDBG_sprint_statuses :
131 #] Utilities :
132 #[ MPI APIs :
133*/
134
135/*
136 * The followings are the MPI APIs with the logging.
137 */
138
139#define MPIDBG_EXTARG const char *file, int line, const char *func
140
141/*
142 #[ MPI_Send :
143*/
144
145static inline int MPIDBG_Send(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
146{
147 MPIDBG_Out("MPI_Send: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
148 int ret = MPI_Send(buf, count, datatype, dest, tag, comm);
149 if ( ret == MPI_SUCCESS ) {
150 MPIDBG_Out("MPI_Send: OK");
151 }
152 else {
153 MPIDBG_Out("MPI_Send: Failed");
154 }
155 return ret;
156}
157#define MPI_Send(...) MPIDBG_Send(__VA_ARGS__, __FILE__, __LINE__, __func__)
158
159/*
160 #] MPI_Send :
161 #[ MPI_Recv :
162*/
163
164static inline int MPIDBG_Recv(void* buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status *status, MPIDBG_EXTARG)
165{
166 MPIDBG_Out("MPI_Recv: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
167 int ret = MPI_Recv(buf, count, datatype, source, tag, comm, status);
168 if ( ret == MPI_SUCCESS ) {
169 int ret_count = 0;
170 MPI_Get_count(status, datatype, &ret_count);
171 MPIDBG_Out("MPI_Recv: OK src=%d dest=%d tag=%d count=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_count);
172 }
173 else {
174 MPIDBG_Out("MPI_Recv: Failed");
175 }
176 return ret;
177}
178#define MPI_Recv(...) MPIDBG_Recv(__VA_ARGS__, __FILE__, __LINE__, __func__)
179
180/*
181 #] MPI_Recv :
182 #[ MPI_Bsend :
183*/
184
185static inline int MPIDBG_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
186{
187 MPIDBG_Out("MPI_Bsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
188 int ret = MPI_Bsend(buf, count, datatype, dest, tag, comm);
189 if ( ret == MPI_SUCCESS ) {
190 MPIDBG_Out("MPI_Bsend: OK");
191 }
192 else {
193 MPIDBG_Out("MPI_Bsend: Failed");
194 }
195 return ret;
196}
197#define MPI_Bsend(...) MPIDBG_Bsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
198
199/*
200 #] MPI_Bsend :
201 #[ MPI_Ssend :
202*/
203
204static inline int MPIDBG_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
205{
206 MPIDBG_Out("MPI_Ssend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
207 int ret = MPI_Ssend(buf, count, datatype, dest, tag, comm);
208 if ( ret == MPI_SUCCESS ) {
209 MPIDBG_Out("MPI_Ssend: OK");
210 }
211 else {
212 MPIDBG_Out("MPI_Ssend: Failed");
213 }
214 return ret;
215}
216#define MPI_Ssend(...) MPIDBG_Ssend(__VA_ARGS__, __FILE__, __LINE__, __func__)
217
218/*
219 #] MPI_Ssend :
220 #[ MPI_Rsend :
221*/
222
223static inline int MPIDBG_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPIDBG_EXTARG)
224{
225 MPIDBG_Out("MPI_Rsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
226 int ret = MPI_Rsend(buf, count, datatype, dest, tag, comm);
227 if ( ret == MPI_SUCCESS ) {
228 MPIDBG_Out("MPI_Rsend: OK");
229 }
230 else {
231 MPIDBG_Out("MPI_Rsend: Failed");
232 }
233 return ret;
234}
235#define MPI_Rsend(...) MPIDBG_Rsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
236
237/*
238 #] MPI_Rsend :
239 #[ MPI_Isend :
240*/
241
242static inline int MPIDBG_Isend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
243{
244 MPIDBG_Out("MPI_Isend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
245 int ret = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
246 if ( ret == MPI_SUCCESS ) {
247 MPIDBG_Out("MPI_Isend: OK");
248 }
249 else {
250 MPIDBG_Out("MPI_Isend: Failed");
251 }
252 return ret;
253}
254#define MPI_Isend(...) MPIDBG_Isend(__VA_ARGS__, __FILE__, __LINE__, __func__)
255
256/*
257 #] MPI_Isend :
258 #[ MPI_Ibsend :
259*/
260
261static inline int MPIDBG_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
262{
263 MPIDBG_Out("MPI_Ibsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
264 int ret = MPI_Ibsend(buf, count, datatype, dest, tag, comm, request);
265 if ( ret == MPI_SUCCESS ) {
266 MPIDBG_Out("MPI_Ibsend: OK");
267 }
268 else {
269 MPIDBG_Out("MPI_Ibsend: Failed");
270 }
271 return ret;
272}
273#define MPI_Ibsend(...) MPIDBG_Ibsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
274
275/*
276 #] MPI_Ibsend :
277 #[ MPI_Issend :
278*/
279
280static inline int MPIDBG_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
281{
282 MPIDBG_Out("MPI_Issend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
283 int ret = MPI_Issend(buf, count, datatype, dest, tag, comm, request);
284 if ( ret == MPI_SUCCESS ) {
285 MPIDBG_Out("MPI_Issend: OK");
286 }
287 else {
288 MPIDBG_Out("MPI_Issend: Failed");
289 }
290 return ret;
291}
292#define MPI_Issend(...) MPIDBG_Issend(__VA_ARGS__, __FILE__, __LINE__, __func__)
293
294/*
295 #] MPI_Issend :
296 #[ MPI_Irsend :
297*/
298
299static inline int MPIDBG_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
300{
301 MPIDBG_Out("MPI_Irsend: src=%d dest=%d tag=%d count=%d", MPIDBG_RANK, dest, tag, count);
302 int ret = MPI_Irsend(buf, count, datatype, dest, tag, comm, request);
303 if ( ret == MPI_SUCCESS ) {
304 MPIDBG_Out("MPI_Irsend: OK");
305 }
306 else {
307 MPIDBG_Out("MPI_Irsend: Failed");
308 }
309 return ret;
310}
311#define MPI_Irsend(...) MPIDBG_Irsend(__VA_ARGS__, __FILE__, __LINE__, __func__)
312
313/*
314 #] MPI_Irsend :
315 #[ MPI_Irecv :
316*/
317
318static inline int MPIDBG_Irecv(void* buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Request *request, MPIDBG_EXTARG)
319{
320 MPIDBG_Out("MPI_Irecv: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
321 int ret = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
322 if ( ret == MPI_SUCCESS ) {
323 MPIDBG_Out("MPI_Irecv: OK dest=%d", MPIDBG_RANK);
324 }
325 else {
326 MPIDBG_Out("MPI_Irecv: Failed");
327 }
328 return ret;
329}
330#define MPI_Irecv(...) MPIDBG_Irecv(__VA_ARGS__, __FILE__, __LINE__, __func__)
331
332/*
333 #] MPI_Irecv :
334 #[ MPI_Wait :
335*/
336
337static inline int MPIDBG_Wait(MPI_Request *request, MPI_Status *status, MPIDBG_EXTARG)
338{
339 char buf[256 * 1]; /* Enough. */
340 MPI_Request old_request = *request;
341 MPIDBG_sprint_requests(buf, 1, request);
342 MPIDBG_Out("MPI_Wait: rank=%d request=%s", MPIDBG_RANK, buf);
343 int ret = MPI_Wait(request, status);
344 if ( ret == MPI_SUCCESS ) {
345 MPIDBG_sprint_statuses(buf, 1, request, &old_request, status);
346 MPIDBG_Out("MPI_Wait: OK rank=%d result=%s", MPIDBG_RANK, buf);
347 }
348 else {
349 MPIDBG_Out("MPI_Wait: Failed");
350 }
351 return ret;
352}
353#define MPI_Wait(...) MPIDBG_Wait(__VA_ARGS__, __FILE__, __LINE__, __func__)
354
355/*
356 #] MPI_Wait :
357 #[ MPI_Test :
358*/
359
360static inline int MPIDBG_Test(MPI_Request *request, int *flag, MPI_Status *status, MPIDBG_EXTARG)
361{
362 char buf[256 * 1]; /* Enough. */
363 MPI_Request old_request = *request;
364 MPIDBG_sprint_requests(buf, 1, request);
365 MPIDBG_Out("MPI_Test: rank=%d request=%s", MPIDBG_RANK, buf);
366 int ret = MPI_Test(request, flag, status);
367 if ( ret == MPI_SUCCESS ) {
368 if ( *flag ) {
369 MPIDBG_sprint_statuses(buf, 1, request, &old_request, status);
370 MPIDBG_Out("MPI_Test: OK rank=%d result=%s", MPIDBG_RANK, buf);
371 }
372 else {
373 MPIDBG_Out("MPI_Test: OK flag=false");
374 }
375 }
376 else {
377 MPIDBG_Out("MPI_Test: Failed");
378 }
379 return ret;
380}
381#define MPI_Test(...) MPIDBG_Test(__VA_ARGS__, __FILE__, __LINE__, __func__)
382
383/*
384 #] MPI_Test :
385 #[ MPI_Waitany :
386*/
387
388static inline int MPIDBG_Waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status, MPIDBG_EXTARG)
389{
390 char buf[256]; /* Enough. */
391 MPI_Request old_requests[count];
392 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
393 MPIDBG_sprint_requests(buf, count, array_of_requests);
394 MPIDBG_Out("MPI_Waitany: rank=%d request=%s", MPIDBG_RANK, buf);
395 int ret = MPI_Waitany(count, array_of_requests, index, status);
396 if ( ret == MPI_SUCCESS ) {
397 MPI_Status statuses[count];
398 statuses[*index] = *status;
399 MPIDBG_sprint_statuses(buf, count, old_requests, array_of_requests, statuses);
400 MPIDBG_Out("MPI_Waitany: OK rank=%d result=%s", MPIDBG_RANK, buf);
401 }
402 else {
403 MPIDBG_Out("MPI_Waitany: Failed");
404 }
405 return ret;
406}
407#define MPI_Waitany(...) MPIDBG_Waitany(__VA_ARGS__, __FILE__, __LINE__, __func__)
408
409/*
410 #] MPI_Waitany :
411 #[ MPI_Testany :
412*/
413
414static inline int MPIDBG_Testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status, MPIDBG_EXTARG)
415{
416 char buf[256]; /* Enough. */
417 MPI_Request old_requests[count];
418 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
419 MPIDBG_sprint_requests(buf, count, array_of_requests);
420 MPIDBG_Out("MPI_Testany: rank=%d request=%s", MPIDBG_RANK, buf);
421 int ret = MPI_Testany(count, array_of_requests, index, flag, status);
422 if ( ret == MPI_SUCCESS ) {
423 if ( *flag ) {
424 MPI_Status statuses[count];
425 statuses[*index] = *status;
426 MPIDBG_sprint_statuses(buf, count, old_requests, array_of_requests, statuses);
427 MPIDBG_Out("MPI_Testany: OK rank=%d result=%s", MPIDBG_RANK, buf);
428 }
429 else {
430 MPIDBG_Out("MPI_Testany: OK flag=false");
431 }
432 }
433 else {
434 MPIDBG_Out("MPI_Testany: Failed");
435 }
436 return ret;
437}
438#define MPI_Testany(...) MPIDBG_Testany(__VA_ARGS__, __FILE__, __LINE__, __func__)
439
440/*
441 #] MPI_Testany :
442 #[ MPI_Waitall :
443*/
444
445static inline int MPIDBG_Waitall(int count, MPI_Request *array_of_requests, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
446{
447 char buf[256 * count]; /* Enough. */
448 MPI_Request old_requests[count];
449 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
450 MPIDBG_sprint_requests(buf, count, array_of_requests);
451 MPIDBG_Out("MPI_Waitall: rank=%d request=%s", MPIDBG_RANK, buf);
452 int ret = MPI_Waitall(count, array_of_requests, array_of_statuses);
453 if ( ret == MPI_SUCCESS ) {
454 MPIDBG_sprint_statuses(buf, count, old_requests, array_of_requests, array_of_statuses);
455 MPIDBG_Out("MPI_Waitall: OK rank=%d result=%s", MPIDBG_RANK, buf);
456 }
457 else {
458 MPIDBG_Out("MPI_Waitall: Failed");
459 }
460 return ret;
461}
462#define MPI_Waitall(...) MPIDBG_Waitall(__VA_ARGS__, __FILE__, __LINE__, __func__)
463
464/*
465 #] MPI_Waitall :
466 #[ MPI_Testall :
467*/
468
469static inline int MPIDBG_Testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
470{
471 char buf[256 * count]; /* Enough. */
472 MPI_Request old_requests[count];
473 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * count);
474 MPIDBG_sprint_requests(buf, count, array_of_requests);
475 MPIDBG_Out("MPI_Testall: rank=%d request=%s", MPIDBG_RANK, buf);
476 int ret = MPI_Testall(count, array_of_requests, flag, array_of_statuses);
477 if ( ret == MPI_SUCCESS ) {
478 if ( *flag ) {
479 MPIDBG_sprint_statuses(buf, count, old_requests, array_of_requests, array_of_statuses);
480 MPIDBG_Out("MPI_Testall: OK rank=%d result=%s", MPIDBG_RANK, buf);
481 }
482 else {
483 MPIDBG_Out("MPI_Testall: OK flag=false");
484 }
485 }
486 else {
487 MPIDBG_Out("MPI_Testall: Failed");
488 }
489 return ret;
490}
491#define MPI_Testall(...) MPIDBG_Testall(__VA_ARGS__, __FILE__, __LINE__, __func__)
492
493/*
494 #] MPI_Testall :
495 #[ MPI_Waitsome :
496*/
497
498static inline int MPIDBG_Waitsome(int incount, MPI_Request *array_of_requests, int *outcount, int *array_of_indices, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
499{
500 char buf[256 * incount]; /* Enough. */
501 MPI_Request old_requests[incount];
502 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * incount);
503 MPIDBG_sprint_requests(buf, incount, array_of_requests);
504 MPIDBG_Out("MPI_Waitsome: rank=%d request=%s", MPIDBG_RANK, buf);
505 int ret = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
506 if ( ret == MPI_SUCCESS ) {
507 MPIDBG_sprint_statuses(buf, incount, old_requests, array_of_requests, array_of_statuses);
508 MPIDBG_Out("MPI_Waitsome: OK rank=%d result=%s", MPIDBG_RANK, buf);
509 }
510 else {
511 MPIDBG_Out("MPI_Waitsome: Failed");
512 }
513 return ret;
514}
515#define MPI_Waitsome(...) MPIDBG_Waitsome(__VA_ARGS__, __FILE__, __LINE__, __func__)
516
517/*
518 #] MPI_Waitsome :
519 #[ MPI_Testsome :
520*/
521
522static inline int MPIDBG_Testsome(int incount, MPI_Request *array_of_requests, int *outcount, int *array_of_indices, MPI_Status *array_of_statuses, MPIDBG_EXTARG)
523{
524 char buf[256 * incount]; /* Enough. */
525 MPI_Request old_requests[incount];
526 memcpy(old_requests, array_of_requests, sizeof(MPI_Request) * incount);
527 MPIDBG_sprint_requests(buf, incount, array_of_requests);
528 MPIDBG_Out("MPI_Testsome: rank=%d request=%s", MPIDBG_RANK, buf);
529 int ret = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
530 if ( ret == MPI_SUCCESS ) {
531 MPIDBG_sprint_statuses(buf, incount, old_requests, array_of_requests, array_of_statuses);
532 MPIDBG_Out("MPI_Testsome: OK rank=%d result=%s", MPIDBG_RANK, buf);
533 }
534 else {
535 MPIDBG_Out("MPI_Testsome: Failed");
536 }
537 return ret;
538}
539#define MPI_Testsome(...) MPIDBG_Testsome(__VA_ARGS__, __FILE__, __LINE__, __func__)
540
541/*
542 #] MPI_Testsome :
543 #[ MPI_Iprobe :
544*/
545
546static inline int MPIDBG_Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status, MPIDBG_EXTARG)
547{
548 MPIDBG_Out("MPI_Iprobe: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
549 int ret = MPI_Iprobe(source, tag, comm, flag, status);
550 if ( ret == MPI_SUCCESS ) {
551 if ( *flag ) {
552 int ret_size = 0;
553 MPI_Get_count(status, MPI_BYTE, &ret_size);
554 MPIDBG_Out("MPI_Iprobe: OK src=%d dest=%d tag=%d size=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_size);
555 }
556 else {
557 MPIDBG_Out("MPI_Iprobe: OK flag=false");
558 }
559 }
560 else {
561 MPIDBG_Out("MPI_Iprobe: Failed");
562 }
563 return ret;
564}
565#define MPI_Iprobe(...) MPIDBG_Iprobe(__VA_ARGS__, __FILE__, __LINE__, __func__)
566
567/*
568 #] MPI_Iprobe :
569 #[ MPI_Probe :
570*/
571
572static inline int MPIDBG_Probe(int source, int tag, MPI_Comm comm, MPI_Status *status, MPIDBG_EXTARG)
573{
574 MPIDBG_Out("MPI_Probe: src=%d dest=%d tag=%d", source, MPIDBG_RANK, tag);
575 int ret = MPI_Probe(source, tag, comm, status);
576 if ( ret == MPI_SUCCESS ) {
577 int ret_size = 0;
578 MPI_Get_count(status, MPI_BYTE, &ret_size);
579 MPIDBG_Out("MPI_Probe: OK src=%d dest=%d tag=%d size=%d", status->MPI_SOURCE, MPIDBG_RANK, status->MPI_TAG, ret_size);
580 }
581 else {
582 MPIDBG_Out("MPI_Probe: Failed");
583 }
584 return ret;
585}
586#define MPI_Probe(...) MPIDBG_Probe(__VA_ARGS__, __FILE__, __LINE__, __func__)
587
588/*
589 #] MPI_Probe :
590 #[ MPI_Cancel :
591*/
592
593static inline int MPIDBG_Cancel(MPI_Request *request, MPIDBG_EXTARG)
594{
595 MPIDBG_Out("MPI_Cancel: rank=%d", MPIDBG_RANK);
596 int ret = MPI_Cancel(request);
597 if ( ret == MPI_SUCCESS ) {
598 MPIDBG_Out("MPI_Cancel: OK");
599 }
600 else {
601 MPIDBG_Out("MPI_Cancel: Failed");
602 }
603 return ret;
604}
605#define MPI_Cancel(...) MPIDBG_Cancel(__VA_ARGS__, __FILE__, __LINE__, __func__)
606
607/*
608 #] MPI_Cancel :
609 #[ MPI_Test_cancelled :
610*/
611
612static inline int MPIDBG_Test_cancelled(MPI_Status *status, int *flag, MPIDBG_EXTARG)
613{
614 MPIDBG_Out("MPI_Test_cancelled: rank=%d", MPIDBG_RANK);
615 int ret = MPI_Test_cancelled(status, flag);
616 if ( ret == MPI_SUCCESS ) {
617 if ( *flag ) {
618 MPIDBG_Out("MPI_Test_cancelled: OK flag=true");
619 }
620 else {
621 MPIDBG_Out("MPI_Test_cancelled: OK flag=false");
622 }
623 }
624 else {
625 MPIDBG_Out("MPI_Test_cancelled: Failed");
626 }
627 return ret;
628}
629#define MPI_Test_cancelled(...) MPIDBG_Test_cancelled(__VA_ARGS__, __FILE__, __LINE__, __func__)
630
631/*
632 #] MPI_Test_cancelled :
633 #[ MPI_Barrier :
634*/
635
636static inline int MPIDBG_Barrier(MPI_Comm comm, MPIDBG_EXTARG)
637{
638 MPIDBG_Out("MPI_Barrier: rank=%d", MPIDBG_RANK);
639 int ret = MPI_Barrier(comm);
640 if ( ret == MPI_SUCCESS ) {
641 MPIDBG_Out("MPI_Barrier: OK");
642 }
643 else {
644 MPIDBG_Out("MPI_Barrier: Failed");
645 }
646 return ret;
647}
648#define MPI_Barrier(...) MPIDBG_Barrier(__VA_ARGS__, __FILE__, __LINE__, __func__)
649
650/*
651 #] MPI_Barrier :
652 #[ MPI_Bcast :
653*/
654
655static inline int MPIDBG_Bcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm, MPIDBG_EXTARG)
656{
657 MPIDBG_Out("MPI_Bcast: root=%d count=%d", root, count);
658 int ret = MPI_Bcast(buffer, count, datatype, root, comm);
659 if ( ret == MPI_SUCCESS ) {
660 MPIDBG_Out("MPI_Bcast: OK");
661 }
662 else {
663 MPIDBG_Out("MPI_Bcast: Failed");
664 }
665 return ret;
666}
667#define MPI_Bcast(...) MPIDBG_Bcast(__VA_ARGS__, __FILE__, __LINE__, __func__)
668
669/*
670 #] MPI_Bcast :
671 #] MPI APIs :
672*/
673
674#endif /* MPIDBG_H */