Arcane  v3.16.0.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
Otf2MpiProfiling.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2025 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 de l'interface IMpiProfiling permettant l'instrumentation */
11/* au format OTF2 . */
12/*---------------------------------------------------------------------------*/
13
14#include "arccore/message_passing_mpi/internal/MessagePassingMpiEnum.h"
15#include "arcane/utils/FatalErrorException.h"
16#include "arcane/std/internal/Otf2MpiProfiling.h"
17
18/*---------------------------------------------------------------------------*/
19/*---------------------------------------------------------------------------*/
20
21namespace Arcane
22{
23using namespace MessagePassing::Mpi;
24
25/*---------------------------------------------------------------------------*/
26/*---------------------------------------------------------------------------*/
27// Helper sur la fonction de taille de type MPI...
28namespace
29{
30using ReturnType = Otf2MpiProfiling::ReturnType;
31
32uint64_t
33getSizeOfMpiType(MPI_Datatype datatype)
34{
35 int tmp_size(0);
36 MPI_Type_size(datatype, &tmp_size);
37 return static_cast<uint64_t>(tmp_size);
38}
39
40}
41
42/*---------------------------------------------------------------------------*/
43/*---------------------------------------------------------------------------*/
45// Note : On ne passe pas le IMpiProfiling a decorer, on appel directement MPI
46// pour eviter un un autre appel inutile. Mais on peut facilement le changer le
47// cas echeant en rajoutant un IMpiProfiling dans les parametres du ctor.
50: m_otf2_wrapper(otf2_wrapper)
51{
52}
53
54/*---------------------------------------------------------------------------*/
55/*---------------------------------------------------------------------------*/
56
59broadcast(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
60{
61 // calcul des tailles des donnees echangees
62 const uint64_t byte_send(m_otf2_wrapper->getMpiRank()==root?getSizeOfMpiType(datatype)*count:0);
63 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank()==root?0:getSizeOfMpiType(datatype)*count);
64
65 _doEventEnter(eMpiName::Bcast);
66
67 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
68
69 int r = MPI_Bcast(buffer, count, datatype, root, comm);
70
71 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
72 OTF2_COLLECTIVE_OP_BCAST, 0 /* comm region */, root,
73 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
74
75 _doEventLeave(eMpiName::Bcast);
76 return _ret(r);
77}
78
79/*---------------------------------------------------------------------------*/
80/*---------------------------------------------------------------------------*/
81
84gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
85 int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
86{
87 // calcul des tailles des donnees echangees
88 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
89 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank()==root?getSizeOfMpiType(recvtype)*recvcount:0);
90
91 _doEventEnter(eMpiName::Gather);
92
93 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
94
95 int r = MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
96
97 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
98 OTF2_COLLECTIVE_OP_GATHER, 0 /* comm region */, root,
99 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
100
101 _doEventLeave(eMpiName::Gather);
102 return _ret(r);
103}
104
105/*---------------------------------------------------------------------------*/
106/*---------------------------------------------------------------------------*/
107
110gatherVariable(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
111 const int *recvcounts, const int *displs, MPI_Datatype recvtype, int root, MPI_Comm comm)
112{
113 // calcul des tailles des donnees echangees
114 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
115 uint64_t byte_recv(0);
116
117 auto f_byte_recv = [&](){for(int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
118 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype);};
119 if (m_otf2_wrapper->getMpiRank() == root)
120 f_byte_recv();
121
122 _doEventEnter(eMpiName::Gatherv);
123
124 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
125
126 int r = MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
127
128 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
129 OTF2_COLLECTIVE_OP_GATHERV, 0 /* comm region */, root,
130 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
131
132 _doEventLeave(eMpiName::Gatherv);
133 return _ret(r);
134}
135
136/*---------------------------------------------------------------------------*/
137/*---------------------------------------------------------------------------*/
138
141allGather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
142 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
143{
144 // calcul des tailles des donnees echangees
145 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
146 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount * m_otf2_wrapper->getMpiNbRank());
147
148 _doEventEnter(eMpiName::Allgather);
149
150 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
151
152 int r = MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
153
154 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
155 OTF2_COLLECTIVE_OP_ALLGATHER, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
156 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
157
158 _doEventLeave(eMpiName::Allgather);
159 return _ret(r);
160}
161
162/*---------------------------------------------------------------------------*/
163/*---------------------------------------------------------------------------*/
164
167allGatherVariable(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
168 const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm)
169{
170 // calcul des tailles des donnees echangees
171 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
172 uint64_t byte_recv(0);
173
174 for (int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
175 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype);
176
177 _doEventEnter(eMpiName::Allgatherv);
178
179 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
180
181 int r = MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
182
183 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
184 OTF2_COLLECTIVE_OP_ALLGATHERV, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
185 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
186
187 _doEventLeave(eMpiName::Allgatherv);
188 return _ret(r);
189}
190
191/*---------------------------------------------------------------------------*/
192/*---------------------------------------------------------------------------*/
193
196scatterVariable(const void *sendbuf, const int *sendcounts, const int *displs,
197 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype,
198 int root, MPI_Comm comm)
199{
200 // calcul des tailles des donnees echangees
201 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount);
202 uint64_t byte_send(0);
203
204 auto f_byte_send = [&](){for(int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
205 byte_send += sendcounts[i];
206 byte_send *= getSizeOfMpiType(recvtype);};
207 if (m_otf2_wrapper->getMpiRank() == root)
208 f_byte_send();
209
210 _doEventEnter(eMpiName::Scatterv);
211
212 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
213
214 int r = MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
215
216 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
217 OTF2_COLLECTIVE_OP_SCATTERV, 0 /* comm region */, root,
218 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
219
220 _doEventLeave(eMpiName::Scatterv);
221 return _ret(r);
222}
223
224/*---------------------------------------------------------------------------*/
225/*---------------------------------------------------------------------------*/
226
229allToAll(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
230 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
231{
232 _doEventEnter(eMpiName::Alltoall);
233
234 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
235
236 int r = MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
237
238 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
239 OTF2_COLLECTIVE_OP_ALLTOALL, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
240 0 /* bytes provided */, 0 /* bytes obtained */);
241
242 _doEventLeave(eMpiName::Alltoall);
243 return _ret(r);
244}
245
246/*---------------------------------------------------------------------------*/
247/*---------------------------------------------------------------------------*/
248
251allToAllVariable(const void *sendbuf, const int *sendcounts, const int *sdispls,
252 MPI_Datatype sendtype, void *recvbuf, const int *recvcounts,
253 const int *rdispls, MPI_Datatype recvtype, MPI_Comm comm)
254{
255 _doEventEnter(eMpiName::Alltoallv);
256
257 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
258
259 int r = MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf, recvcounts, rdispls, recvtype, comm);
260
261 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
262 OTF2_COLLECTIVE_OP_ALLTOALLV, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
263 0 /* bytes provided */, 0 /* bytes obtained */);
264
265 _doEventLeave(eMpiName::Alltoallv);
266 return _ret(r);
267}
268
269/*---------------------------------------------------------------------------*/
270/*---------------------------------------------------------------------------*/
271
274barrier(MPI_Comm comm)
275{
276 _doEventEnter(eMpiName::Barrier);
277
278 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
279
280 int r = MPI_Barrier(comm);
281
282 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
283 OTF2_COLLECTIVE_OP_BARRIER, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
284 0 /* bytes provided */, 0 /* bytes obtained */);
285
286 _doEventLeave(eMpiName::Barrier);
287 return _ret(r);
288}
289
290/*---------------------------------------------------------------------------*/
291/*---------------------------------------------------------------------------*/
292
295reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
296 MPI_Op op, int root, MPI_Comm comm)
297{
298 _doEventEnter(eMpiName::Reduce);
299
300 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
301
302 int r = MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
303
304 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
305 OTF2_COLLECTIVE_OP_REDUCE, 0 /* comm region */, root,
306 0 /* bytes provided */, 0 /* bytes obtained */);
307
308 _doEventLeave(eMpiName::Reduce);
309 return _ret(r);
310}
311
312/*---------------------------------------------------------------------------*/
313/*---------------------------------------------------------------------------*/
314
317allReduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
318{
319 // calcul des tailles des donnees echangees
320 const uint64_t byte_send(getSizeOfMpiType(datatype) * count);
321 const uint64_t byte_recv(getSizeOfMpiType(datatype) * count);
322
323 _doEventEnter(eMpiName::Allreduce);
324
325 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
326
327 int r = MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
328
329 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
330 OTF2_COLLECTIVE_OP_ALLREDUCE, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
331 byte_send /* bytes provided */, byte_recv /* bytes obtained */);
332
333 _doEventLeave(eMpiName::Allreduce);
334 return _ret(r);
335}
336
337/*---------------------------------------------------------------------------*/
338/*---------------------------------------------------------------------------*/
339
342scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
343{
344 _doEventEnter(eMpiName::Scan);
345
346 OTF2_EvtWriter_MpiCollectiveBegin(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime());
347
348 int r = MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
349
350 int type_size;
351 MPI_Type_size(datatype, &type_size);
352 const uint64_t size = static_cast<uint64_t>(type_size) * static_cast<uint64_t>(count);
353 OTF2_EvtWriter_MpiCollectiveEnd(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
354 OTF2_COLLECTIVE_OP_SCAN, 0 /* comm region */, OTF2_UNDEFINED_UINT32 /* root */,
355 size /* bytes provided */, size /* bytes obtained */);
356
357 _doEventLeave(eMpiName::Scan);
358 return _ret(r);
359}
360
361/*---------------------------------------------------------------------------*/
362/*---------------------------------------------------------------------------*/
363
366sendRecv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest,
367 int sendtag, void *recvbuf, int recvcount, MPI_Datatype recvtype,
368 int source, int recvtag, MPI_Comm comm, MPI_Status *status)
369{
370 _doEventEnter(eMpiName::Sendrecv);
371
372 int r = MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype,
373 source, recvtag, comm, status);
374
375 _doEventLeave(eMpiName::Sendrecv);
376 return _ret(r);
377}
378
379/*---------------------------------------------------------------------------*/
380/*---------------------------------------------------------------------------*/
381
384iSend(const void *buf, int count, MPI_Datatype datatype, int dest, int tag,
385 MPI_Comm comm, MPI_Request *request)
386{
387 _doEventEnter(eMpiName::Isend);
388
389 int r = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
390 // Comme OTF2 a besoin d'un entier, utilise la valeur fortran de la requête qui
391 // est garantie être un entier (contrairement à MPI_Request)
392 uint64_t request_id = MPI_Request_c2f(*request);
393 OTF2_EvtWriter_MpiIsend(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
394 dest, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count, request_id);
395
396 _doEventLeave(eMpiName::Isend);
397 return _ret(r);
398}
399
400/*---------------------------------------------------------------------------*/
401/*---------------------------------------------------------------------------*/
402
405send(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm)
406{
407 _doEventEnter(eMpiName::Send);
408
409 int r = MPI_Send(buf, count, datatype, dest, tag, comm);
410
411 OTF2_EvtWriter_MpiSend(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
412 dest, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count);
413
414 _doEventLeave(eMpiName::Send);
415 return _ret(r);
416}
417
418/*---------------------------------------------------------------------------*/
419/*---------------------------------------------------------------------------*/
420
423iRecv(void *buf, int count, MPI_Datatype datatype, int source, int tag,
424 MPI_Comm comm, MPI_Request *request)
425{
426 _doEventEnter(eMpiName::Irecv);
427
428 int r = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
429
430 // TODO: arriver a faire fonctionner ce truc...
431// OTF2_EvtWriter_MpiIrecv(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
432// source, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count, *request);
433 _doEventLeave(eMpiName::Irecv);
434 return _ret(r);
435}
436
437/*---------------------------------------------------------------------------*/
438/*---------------------------------------------------------------------------*/
439
442recv(void *buf, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status *status)
443{
444 _doEventEnter(eMpiName::Recv);
445 int r = MPI_Recv(buf, count, datatype, source, tag, comm, status);
446
447 OTF2_EvtWriter_MpiRecv(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
448 source, 0 /* comm region */, tag, getSizeOfMpiType(datatype) * count);
449
450 _doEventLeave(eMpiName::Recv);
451 return _ret(r);
452}
453
454/*---------------------------------------------------------------------------*/
455/*---------------------------------------------------------------------------*/
456
459test(MPI_Request *request, int *flag, MPI_Status *status)
460{
461 _doEventEnter(eMpiName::Test);
462 if (!MPI_Test(request, flag, status)){
463 uint64_t request_id = MPI_Request_c2f(*request);
464 OTF2_EvtWriter_MpiRequestTest(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(), request_id);
465 }
466 _doEventLeave(eMpiName::Test);
467 return _ret(0);
468}
469
470/*---------------------------------------------------------------------------*/
471/*---------------------------------------------------------------------------*/
472
475probe(int source, int tag, MPI_Comm comm, MPI_Status *status)
476{
477 _doEventEnter(eMpiName::Probe);
478 int r = MPI_Probe(source, tag, comm, status);
479 _doEventLeave(eMpiName::Probe);
480 return _ret(r);
481}
482
483/*---------------------------------------------------------------------------*/
484/*---------------------------------------------------------------------------*/
485
488getCount(const MPI_Status *status, MPI_Datatype datatype, int *count)
489{
490 _doEventEnter(eMpiName::Get_count);
491 int r = MPI_Get_count(status, datatype, count);
492 _doEventLeave(eMpiName::Get_count);
493 return _ret(r);
494}
495
496/*---------------------------------------------------------------------------*/
497/*---------------------------------------------------------------------------*/
498
501wait(MPI_Request *request, MPI_Status *status)
502{
503 _doEventEnter(eMpiName::Wait);
504 int r = MPI_Wait(request, status);
505 _doEventLeave(eMpiName::Wait);
506 return _ret(r);
507}
508
509/*---------------------------------------------------------------------------*/
510/*---------------------------------------------------------------------------*/
511
514waitAll(int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[])
515{
516 _doEventEnter(eMpiName::Waitall);
517 int r = MPI_Waitall(count, array_of_requests, array_of_statuses);
518 _doEventLeave(eMpiName::Waitall);
519 return _ret(r);
520}
521
522/*---------------------------------------------------------------------------*/
523/*---------------------------------------------------------------------------*/
524
527testSome(int incount, MPI_Request array_of_requests[], int *outcount,
528 int array_of_indices[], MPI_Status array_of_statuses[])
529{
530 _doEventEnter(eMpiName::Testsome);
531 int r = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
532 _doEventLeave(eMpiName::Testsome);
533 return _ret(r);
534}
535
536/*---------------------------------------------------------------------------*/
537/*---------------------------------------------------------------------------*/
538
541waitSome(int incount, MPI_Request array_of_requests[], int *outcount,
542 int array_of_indices[], MPI_Status array_of_statuses[])
543{
544 _doEventEnter(eMpiName::Waitsome);
545 int r = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
546 _doEventLeave(eMpiName::Waitsome);
547 return _ret(r);
548}
549
550/*---------------------------------------------------------------------------*/
551/*---------------------------------------------------------------------------*/
552
553void Otf2MpiProfiling::
554_doEventEnter(eMpiName event_name)
555{
556 OTF2_EvtWriter_Enter(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
557 static_cast<uint32_t>(event_name));
558}
559
560/*---------------------------------------------------------------------------*/
561/*---------------------------------------------------------------------------*/
562
563void Otf2MpiProfiling::
564_doEventLeave(eMpiName event_name)
565{
566 OTF2_EvtWriter_Leave(m_otf2_wrapper->getEventWriter(), NULL, Otf2LibWrapper::getTime(),
567 static_cast<uint32_t>(event_name));
568
569}
570
571/*---------------------------------------------------------------------------*/
572/*---------------------------------------------------------------------------*/
573
574} // namespace Arcane
575
576/*---------------------------------------------------------------------------*/
577/*---------------------------------------------------------------------------*/
Classe d'encapsulation des fonctions de la librairie OTF2.
OTF2_EvtWriter * getEventWriter()
static OTF2_TimeStamp getTime()
Méthode interne statique pour recuperer le 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)
Constructeur.
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 -*-