Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
Otf2MpiProfiling.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* Otf2MpiProfiling.cc (C) 2000-2025 */
9/* */
10/* Implementation of the IMpiProfiling interface allowing instrumentation */
11/* in OTF2 format . */
12/*---------------------------------------------------------------------------*/
13/*---------------------------------------------------------------------------*/
14
15#include "arccore/message_passing_mpi/internal/MessagePassingMpiEnum.h"
16#include "arcane/utils/FatalErrorException.h"
17#include "arcane/std/internal/Otf2MpiProfiling.h"
18
19/*---------------------------------------------------------------------------*/
20/*---------------------------------------------------------------------------*/
21
22namespace Arcane
23{
24using namespace MessagePassing::Mpi;
25
26/*---------------------------------------------------------------------------*/
27/*---------------------------------------------------------------------------*/
28
29// Helper for the MPI type size function...
30namespace
31{
32 using ReturnType = Otf2MpiProfiling::ReturnType;
33
34 uint64_t
35 getSizeOfMpiType(MPI_Datatype datatype)
36 {
37 int tmp_size(0);
38 MPI_Type_size(datatype, &tmp_size);
39 return static_cast<uint64_t>(tmp_size);
40 }
41
42} // namespace
43
44/*---------------------------------------------------------------------------*/
45/*---------------------------------------------------------------------------*/
46
48// Note: We do not pass the IMpiProfiling to be decorated; we call MPI directly
49// to avoid another unnecessary call. But we can easily change this if necessary
50// by adding an IMpiProfiling to the ctor parameters.
53: m_otf2_wrapper(otf2_wrapper)
54{
55}
56
57/*---------------------------------------------------------------------------*/
58/*---------------------------------------------------------------------------*/
59
62broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
63{
64 // calculation of the sizes of the exchanged data
65 const uint64_t byte_send(m_otf2_wrapper->getMpiRank() == root ? getSizeOfMpiType(datatype) * count : 0);
66 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank() == root ? 0 : getSizeOfMpiType(datatype) * count);
67
68 _doEventEnter(eMpiName::Bcast);
69
70 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
71
72 int r = MPI_Bcast(buffer, count, datatype, root, comm);
73
74 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
75 OTF2_COLLECTIVE_OP_BCAST, 0 /* comm region */, root,
76 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
77
78 _doEventLeave(eMpiName::Bcast);
79 return _ret(r);
80}
81
82/*---------------------------------------------------------------------------*/
83/*---------------------------------------------------------------------------*/
84
87gather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf,
88 int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
89{
90 // calculation of the sizes of the exchanged data
91 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
92 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank() == root ? getSizeOfMpiType(recvtype) * recvcount : 0);
93
94 _doEventEnter(eMpiName::Gather);
95
96 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
97
98 int r = MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
99
100 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
101 OTF2_COLLECTIVE_OP_GATHER, 0 /* comm region */, root,
102 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
103
104 _doEventLeave(eMpiName::Gather);
105 return _ret(r);
106}
107
108/*---------------------------------------------------------------------------*/
109/*---------------------------------------------------------------------------*/
110
113gatherVariable(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf,
114 const int* recvcounts, const int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm)
115{
116 // calculation of the sizes of the exchanged data
117 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
118 uint64_t byte_recv(0);
119
120 auto f_byte_recv = [&]() {for(int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
121 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype); };
122 if (m_otf2_wrapper->getMpiRank() == root)
123 f_byte_recv();
124
125 _doEventEnter(eMpiName::Gatherv);
126
127 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
128
129 int r = MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
130
131 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
132 OTF2_COLLECTIVE_OP_GATHERV, 0 /* comm region */, root,
133 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
134
135 _doEventLeave(eMpiName::Gatherv);
136 return _ret(r);
137}
138
139/*---------------------------------------------------------------------------*/
140/*---------------------------------------------------------------------------*/
141
144allGather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf,
145 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
146{
147 // calculation of the sizes of the exchanged data
148 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
149 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount * m_otf2_wrapper->getMpiNbRank());
150
151 _doEventEnter(eMpiName::Allgather);
152
153 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
154
155 int r = MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
156
157 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
158 OTF2_COLLECTIVE_OP_ALLGATHER, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
159 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
160
161 _doEventLeave(eMpiName::Allgather);
162 return _ret(r);
163}
164
165/*---------------------------------------------------------------------------*/
166/*---------------------------------------------------------------------------*/
167
170allGatherVariable(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf,
171 const int* recvcounts, const int* displs, MPI_Datatype recvtype, MPI_Comm comm)
172{
173 // calculation of the sizes of the exchanged data
174 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
175 uint64_t byte_recv(0);
176
177 for (int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
178 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype);
179
180 _doEventEnter(eMpiName::Allgatherv);
181
182 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
183
184 int r = MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
185
186 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
187 OTF2_COLLECTIVE_OP_ALLGATHERV, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
188 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
189
190 _doEventLeave(eMpiName::Allgatherv);
191 return _ret(r);
192}
193
194/*---------------------------------------------------------------------------*/
195/*---------------------------------------------------------------------------*/
196
199scatterVariable(const void* sendbuf, const int* sendcounts, const int* displs,
200 MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype,
201 int root, MPI_Comm comm)
202{
203 // calculation of the sizes of the exchanged data
204 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount);
205 uint64_t byte_send(0);
206
207 auto f_byte_send = [&]() {for(int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
208 byte_send += sendcounts[i];
209 byte_send *= getSizeOfMpiType(recvtype); };
210 if (m_otf2_wrapper->getMpiRank() == root)
211 f_byte_send();
212
213 _doEventEnter(eMpiName::Scatterv);
214
215 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
216
217 int r = MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
218
219 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
220 OTF2_COLLECTIVE_OP_SCATTERV, 0 /* comm region */, root,
221 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
222
223 _doEventLeave(eMpiName::Scatterv);
224 return _ret(r);
225}
226
227/*---------------------------------------------------------------------------*/
228/*---------------------------------------------------------------------------*/
229
232allToAll(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf,
233 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
234{
235 _doEventEnter(eMpiName::Alltoall);
236
237 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
238
239 int r = MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
240
241 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
242 OTF2_COLLECTIVE_OP_ALLTOALL, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
243 0 /* bytes provided */, 0 /* bytes obtained */);
244
245 _doEventLeave(eMpiName::Alltoall);
246 return _ret(r);
247}
248
249/*---------------------------------------------------------------------------*/
250/*---------------------------------------------------------------------------*/
251
254allToAllVariable(const void* sendbuf, const int* sendcounts, const int* sdispls,
255 MPI_Datatype sendtype, void* recvbuf, const int* recvcounts,
256 const int* rdispls, MPI_Datatype recvtype, MPI_Comm comm)
257{
258 _doEventEnter(eMpiName::Alltoallv);
259
260 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
261
262 int r = MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf, recvcounts, rdispls, recvtype, comm);
263
264 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
265 OTF2_COLLECTIVE_OP_ALLTOALLV, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
266 0 /* bytes provided */, 0 /* bytes obtained */);
267
268 _doEventLeave(eMpiName::Alltoallv);
269 return _ret(r);
270}
271
272/*---------------------------------------------------------------------------*/
273/*---------------------------------------------------------------------------*/
274
277barrier(MPI_Comm comm)
278{
279 _doEventEnter(eMpiName::Barrier);
280
281 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
282
283 int r = MPI_Barrier(comm);
284
285 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
286 OTF2_COLLECTIVE_OP_BARRIER, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
287 0 /* bytes provided */, 0 /* bytes obtained */);
288
289 _doEventLeave(eMpiName::Barrier);
290 return _ret(r);
291}
292
293/*---------------------------------------------------------------------------*/
294/*---------------------------------------------------------------------------*/
295
298reduce(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype,
299 MPI_Op op, int root, MPI_Comm comm)
300{
301 _doEventEnter(eMpiName::Reduce);
302
303 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
304
305 int r = MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
306
307 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
308 OTF2_COLLECTIVE_OP_REDUCE, 0 /* comm region */, root,
309 0 /* bytes provided */, 0 /* bytes obtained */);
310
311 _doEventLeave(eMpiName::Reduce);
312 return _ret(r);
313}
314
315/*---------------------------------------------------------------------------*/
316/*---------------------------------------------------------------------------*/
317
320allReduce(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
321{
322 // calculation of the sizes of the exchanged data
323 const uint64_t byte_send(getSizeOfMpiType(datatype) * count);
324 const uint64_t byte_recv(getSizeOfMpiType(datatype) * count);
325
326 _doEventEnter(eMpiName::Allreduce);
327
328 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
329
330 int r = MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
331
332 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
333 OTF2_COLLECTIVE_OP_ALLREDUCE, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
334 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
335
336 _doEventLeave(eMpiName::Allreduce);
337 return _ret(r);
338}
339
340/*---------------------------------------------------------------------------*/
341/*---------------------------------------------------------------------------*/
342
345scan(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
346{
347 _doEventEnter(eMpiName::Scan);
348
349 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
350
351 int r = MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
352
353 int type_size;
354 MPI_Type_size(datatype, &type_size);
355 const uint64_t size = static_cast<uint64_t>(type_size) * static_cast<uint64_t>(count);
356 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
357 OTF2_COLLECTIVE_OP_SCAN, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
358 size /* bytes provided */, size /* bytes obtained */);
359
360 _doEventLeave(eMpiName::Scan);
361 return _ret(r);
362}
363
364/*---------------------------------------------------------------------------*/
365/*---------------------------------------------------------------------------*/
366
369sendRecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dest,
370 int sendtag, void* recvbuf, int recvcount, MPI_Datatype recvtype,
371 int source, int recvtag, MPI_Comm comm, MPI_Status* status)
372{
373 _doEventEnter(eMpiName::Sendrecv);
374
375 int r = MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype,
376 source, recvtag, comm, status);
377
378 _doEventLeave(eMpiName::Sendrecv);
379 return _ret(r);
380}
381
382/*---------------------------------------------------------------------------*/
383/*---------------------------------------------------------------------------*/
384
387iSend(const void* buf, int count, MPI_Datatype datatype, int dest, int tag,
388 MPI_Comm comm, MPI_Request* request)
389{
390 _doEventEnter(eMpiName::Isend);
391
392 int r = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
393 // Since OTF2 requires an integer, use the Fortran value of the request which
394 // is guaranteed to be an integer (unlike MPI_Request)
395 uint64_t request_id = MPI_Request_c2f(*request);
396 OTF2_EvtWriter_MpiIsend(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
397 dest, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count, request_id);
398
399 _doEventLeave(eMpiName::Isend);
400 return _ret(r);
401}
402
403/*---------------------------------------------------------------------------*/
404/*---------------------------------------------------------------------------*/
405
408send(const void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm)
409{
410 _doEventEnter(eMpiName::Send);
411
412 int r = MPI_Send(buf, count, datatype, dest, tag, comm);
413
414 OTF2_EvtWriter_MpiSend(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
415 dest, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count);
416
417 _doEventLeave(eMpiName::Send);
418 return _ret(r);
419}
420
421/*---------------------------------------------------------------------------*/
422/*---------------------------------------------------------------------------*/
423
426iRecv(void* buf, int count, MPI_Datatype datatype, int source, int tag,
427 MPI_Comm comm, MPI_Request* request)
428{
429 _doEventEnter(eMpiName::Irecv);
430
431 int r = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
432
433 // TODO: figure out how to make this thing work...
434 // OTF2_EvtWriter_MpiIrecv(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
435 // source, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count, *request);
436 _doEventLeave(eMpiName::Irecv);
437 return _ret(r);
438}
439
440/*---------------------------------------------------------------------------*/
441/*---------------------------------------------------------------------------*/
442
445recv(void* buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status)
446{
447 _doEventEnter(eMpiName::Recv);
448 int r = MPI_Recv(buf, count, datatype, source, tag, comm, status);
449
450 OTF2_EvtWriter_MpiRecv(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
451 source, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count);
452
453 _doEventLeave(eMpiName::Recv);
454 return _ret(r);
455}
456
457/*---------------------------------------------------------------------------*/
458/*---------------------------------------------------------------------------*/
459
462test(MPI_Request* request, int* flag, MPI_Status* status)
463{
464 _doEventEnter(eMpiName::Test);
465 if (!MPI_Test(request, flag, status)) {
466 uint64_t request_id = MPI_Request_c2f(*request);
467 OTF2_EvtWriter_MpiRequestTest(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(), request_id);
468 }
469 _doEventLeave(eMpiName::Test);
470 return _ret(0);
471}
472
473/*---------------------------------------------------------------------------*/
474/*---------------------------------------------------------------------------*/
475
478probe(int source, int tag, MPI_Comm comm, MPI_Status* status)
479{
480 _doEventEnter(eMpiName::Probe);
481 int r = MPI_Probe(source, tag, comm, status);
482 _doEventLeave(eMpiName::Probe);
483 return _ret(r);
484}
485
486/*---------------------------------------------------------------------------*/
487/*---------------------------------------------------------------------------*/
488
491getCount(const MPI_Status* status, MPI_Datatype datatype, int* count)
492{
493 _doEventEnter(eMpiName::Get_count);
494 int r = MPI_Get_count(status, datatype, count);
495 _doEventLeave(eMpiName::Get_count);
496 return _ret(r);
497}
498
499/*---------------------------------------------------------------------------*/
500/*---------------------------------------------------------------------------*/
501
504wait(MPI_Request* request, MPI_Status* status)
505{
506 _doEventEnter(eMpiName::Wait);
507 int r = MPI_Wait(request, status);
508 _doEventLeave(eMpiName::Wait);
509 return _ret(r);
510}
511
512/*---------------------------------------------------------------------------*/
513/*---------------------------------------------------------------------------*/
514
517waitAll(int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[])
518{
519 _doEventEnter(eMpiName::Waitall);
520 int r = MPI_Waitall(count, array_of_requests, array_of_statuses);
521 _doEventLeave(eMpiName::Waitall);
522 return _ret(r);
523}
524
525/*---------------------------------------------------------------------------*/
526/*---------------------------------------------------------------------------*/
527
530testSome(int incount, MPI_Request array_of_requests[], int* outcount,
531 int array_of_indices[], MPI_Status array_of_statuses[])
532{
533 _doEventEnter(eMpiName::Testsome);
534 int r = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
535 _doEventLeave(eMpiName::Testsome);
536 return _ret(r);
537}
538
539/*---------------------------------------------------------------------------*/
540/*---------------------------------------------------------------------------*/
541
544waitSome(int incount, MPI_Request array_of_requests[], int* outcount,
545 int array_of_indices[], MPI_Status array_of_statuses[])
546{
547 _doEventEnter(eMpiName::Waitsome);
548 int r = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
549 _doEventLeave(eMpiName::Waitsome);
550 return _ret(r);
551}
552
553/*---------------------------------------------------------------------------*/
554/*---------------------------------------------------------------------------*/
555
556void Otf2MpiProfiling::
557_doEventEnter(eMpiName event_name)
558{
559 OTF2_EvtWriter_Enter(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
560 static_cast<uint32_t>(event_name));
561}
562
563/*---------------------------------------------------------------------------*/
564/*---------------------------------------------------------------------------*/
565
566void Otf2MpiProfiling::
567_doEventLeave(eMpiName event_name)
568{
569 OTF2_EvtWriter_Leave(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
570 static_cast<uint32_t>(event_name));
571}
572
573/*---------------------------------------------------------------------------*/
574/*---------------------------------------------------------------------------*/
575
576} // namespace Arcane
577
578/*---------------------------------------------------------------------------*/
579/*---------------------------------------------------------------------------*/
Wrapper class for OTF2 library functions.
OTF2_EvtWriter * getEventWriter()
static OTF2_TimeStamp getTime()
Internal static method to retrieve the timestamp.
ReturnType allGatherVariable(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm) final
MPI_Allgatherv.
ReturnType allGather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) final
MPI_Allgather.
ReturnType waitAll(int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[]) final
MPI_Waitall.
ReturnType testSome(int incount, MPI_Request array_of_requests[], int *outcount, int array_of_indices[], MPI_Status array_of_statuses[]) final
MPI_Testsome.
ReturnType probe(int source, int tag, MPI_Comm comm, MPI_Status *status) final
MPI_Probe.
Otf2MpiProfiling(Otf2LibWrapper *otf2_wrapper)
Constructor.
ReturnType barrier(MPI_Comm comm) final
MPI_Barrier.
ReturnType sendRecv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm, MPI_Status *status) final
MPI_Sendrecv.
ReturnType recv(void *buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status *status) final
MPI_recv.
ReturnType reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) final
MPI_Reduce.
ReturnType scatterVariable(const void *sendbuf, const int *sendcounts, const int *displs, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) final
MPI_Scatterv.
ReturnType allToAllVariable(const void *sendbuf, const int *sendcounts, const int *sdispls, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *rdispls, MPI_Datatype recvtype, MPI_Comm comm) final
MPI_Alltoallv.
ReturnType allToAll(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) final
MPI_Alltoall.
ReturnType gatherVariable(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, int root, MPI_Comm comm) final
MPI_Gatherv.
ReturnType scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) final
MPI_Scan.
ReturnType waitSome(int incount, MPI_Request array_of_requests[], int *outcount, int array_of_indices[], MPI_Status array_of_statuses[]) final
MPI_Waitsome.
ReturnType getCount(const MPI_Status *status, MPI_Datatype datatype, int *count) final
MPI_Get_count.
ReturnType broadcast(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) final
MPI_Bcast.
ReturnType gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) final
MPI_Gather.
ReturnType wait(MPI_Request *request, MPI_Status *status) final
MPI_Wait.
ReturnType test(MPI_Request *request, int *flag, MPI_Status *status) final
MPI_Test.
ReturnType send(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) final
MPI_Send.
ReturnType allReduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) final
MPI_Allreduce.
ReturnType iSend(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request) final
MPI_Isend.
ReturnType iRecv(void *buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Request *request) final
MPI_Irecv.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --