14#include "arccore/message_passing_mpi/internal/MessagePassingMpiEnum.h"
15#include "arcane/utils/FatalErrorException.h"
16#include "arcane/std/internal/Otf2MpiProfiling.h"
23using namespace MessagePassing::Mpi;
30using ReturnType = Otf2MpiProfiling::ReturnType;
33getSizeOfMpiType(MPI_Datatype datatype)
36 MPI_Type_size(datatype, &tmp_size);
37 return static_cast<uint64_t
>(tmp_size);
50: m_otf2_wrapper(otf2_wrapper)
59broadcast(
void *buffer,
int count, MPI_Datatype datatype,
int root, MPI_Comm comm)
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);
65 _doEventEnter(eMpiName::Bcast);
69 int r = MPI_Bcast(buffer, count, datatype, root, comm);
72 OTF2_COLLECTIVE_OP_BCAST, 0 , root,
73 byte_send , byte_recv );
75 _doEventLeave(eMpiName::Bcast);
84gather(
const void *sendbuf,
int sendcount, MPI_Datatype sendtype,
void *recvbuf,
85 int recvcount, MPI_Datatype recvtype,
int root, MPI_Comm comm)
88 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
89 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank()==root?getSizeOfMpiType(recvtype)*recvcount:0);
91 _doEventEnter(eMpiName::Gather);
95 int r = MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
98 OTF2_COLLECTIVE_OP_GATHER, 0 , root,
99 byte_send , byte_recv );
101 _doEventLeave(eMpiName::Gather);
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)
114 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
115 uint64_t byte_recv(0);
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)
122 _doEventEnter(eMpiName::Gatherv);
126 int r = MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
129 OTF2_COLLECTIVE_OP_GATHERV, 0 , root,
130 byte_send , byte_recv );
132 _doEventLeave(eMpiName::Gatherv);
141allGather(
const void *sendbuf,
int sendcount, MPI_Datatype sendtype,
void *recvbuf,
142 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
145 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
146 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount * m_otf2_wrapper->getMpiNbRank());
148 _doEventEnter(eMpiName::Allgather);
152 int r = MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
155 OTF2_COLLECTIVE_OP_ALLGATHER, 0 , OTF2_UNDEFINED_UINT32 ,
156 byte_send , byte_recv );
158 _doEventLeave(eMpiName::Allgather);
167allGatherVariable(
const void *sendbuf,
int sendcount, MPI_Datatype sendtype,
void *recvbuf,
168 const int *recvcounts,
const int *displs, MPI_Datatype recvtype, MPI_Comm comm)
171 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
172 uint64_t byte_recv(0);
174 for (
int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
175 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype);
177 _doEventEnter(eMpiName::Allgatherv);
181 int r = MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
184 OTF2_COLLECTIVE_OP_ALLGATHERV, 0 , OTF2_UNDEFINED_UINT32 ,
185 byte_send , byte_recv );
187 _doEventLeave(eMpiName::Allgatherv);
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)
201 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount);
202 uint64_t byte_send(0);
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)
210 _doEventEnter(eMpiName::Scatterv);
214 int r = MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
217 OTF2_COLLECTIVE_OP_SCATTERV, 0 , root,
218 byte_send , byte_recv );
220 _doEventLeave(eMpiName::Scatterv);
229allToAll(
const void *sendbuf,
int sendcount, MPI_Datatype sendtype,
void *recvbuf,
230 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
232 _doEventEnter(eMpiName::Alltoall);
236 int r = MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
239 OTF2_COLLECTIVE_OP_ALLTOALL, 0 , OTF2_UNDEFINED_UINT32 ,
242 _doEventLeave(eMpiName::Alltoall);
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)
255 _doEventEnter(eMpiName::Alltoallv);
259 int r = MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf, recvcounts, rdispls, recvtype, comm);
262 OTF2_COLLECTIVE_OP_ALLTOALLV, 0 , OTF2_UNDEFINED_UINT32 ,
265 _doEventLeave(eMpiName::Alltoallv);
276 _doEventEnter(eMpiName::Barrier);
280 int r = MPI_Barrier(comm);
283 OTF2_COLLECTIVE_OP_BARRIER, 0 , OTF2_UNDEFINED_UINT32 ,
286 _doEventLeave(eMpiName::Barrier);
295reduce(
const void *sendbuf,
void *recvbuf,
int count, MPI_Datatype datatype,
296 MPI_Op op,
int root, MPI_Comm comm)
298 _doEventEnter(eMpiName::Reduce);
302 int r = MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
305 OTF2_COLLECTIVE_OP_REDUCE, 0 , root,
308 _doEventLeave(eMpiName::Reduce);
317allReduce(
const void *sendbuf,
void *recvbuf,
int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
320 const uint64_t byte_send(getSizeOfMpiType(datatype) * count);
321 const uint64_t byte_recv(getSizeOfMpiType(datatype) * count);
323 _doEventEnter(eMpiName::Allreduce);
327 int r = MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
330 OTF2_COLLECTIVE_OP_ALLREDUCE, 0 , OTF2_UNDEFINED_UINT32 ,
331 byte_send , byte_recv );
333 _doEventLeave(eMpiName::Allreduce);
342scan(
const void *sendbuf,
void *recvbuf,
int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
344 _doEventEnter(eMpiName::Scan);
348 int r = MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
351 MPI_Type_size(datatype, &type_size);
352 const uint64_t size =
static_cast<uint64_t
>(type_size) *
static_cast<uint64_t
>(count);
354 OTF2_COLLECTIVE_OP_SCAN, 0 , OTF2_UNDEFINED_UINT32 ,
357 _doEventLeave(eMpiName::Scan);
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)
370 _doEventEnter(eMpiName::Sendrecv);
372 int r = MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype,
373 source, recvtag, comm, status);
375 _doEventLeave(eMpiName::Sendrecv);
384iSend(
const void *buf,
int count, MPI_Datatype datatype,
int dest,
int tag,
385 MPI_Comm comm, MPI_Request *request)
387 _doEventEnter(eMpiName::Isend);
389 int r = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
392 uint64_t request_id = MPI_Request_c2f(*request);
394 dest, 0 , tag, getSizeOfMpiType(datatype) * count, request_id);
396 _doEventLeave(eMpiName::Isend);
405send(
const void *buf,
int count, MPI_Datatype datatype,
int dest,
int tag, MPI_Comm comm)
407 _doEventEnter(eMpiName::Send);
409 int r = MPI_Send(buf, count, datatype, dest, tag, comm);
412 dest, 0 , tag, getSizeOfMpiType(datatype) * count);
414 _doEventLeave(eMpiName::Send);
423iRecv(
void *buf,
int count, MPI_Datatype datatype,
int source,
int tag,
424 MPI_Comm comm, MPI_Request *request)
426 _doEventEnter(eMpiName::Irecv);
428 int r = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
433 _doEventLeave(eMpiName::Irecv);
442recv(
void *buf,
int count, MPI_Datatype datatype,
int source,
int tag, MPI_Comm comm, MPI_Status *status)
444 _doEventEnter(eMpiName::Recv);
445 int r = MPI_Recv(buf, count, datatype, source, tag, comm, status);
448 source, 0 , tag, getSizeOfMpiType(datatype) * count);
450 _doEventLeave(eMpiName::Recv);
459test(MPI_Request *request,
int *flag, MPI_Status *status)
461 _doEventEnter(eMpiName::Test);
462 if (!MPI_Test(request, flag, status)){
463 uint64_t request_id = MPI_Request_c2f(*request);
466 _doEventLeave(eMpiName::Test);
475probe(
int source,
int tag, MPI_Comm comm, MPI_Status *status)
477 _doEventEnter(eMpiName::Probe);
478 int r = MPI_Probe(source, tag, comm, status);
479 _doEventLeave(eMpiName::Probe);
488getCount(
const MPI_Status *status, MPI_Datatype datatype,
int *count)
490 _doEventEnter(eMpiName::Get_count);
491 int r = MPI_Get_count(status, datatype, count);
492 _doEventLeave(eMpiName::Get_count);
501wait(MPI_Request *request, MPI_Status *status)
503 _doEventEnter(eMpiName::Wait);
504 int r = MPI_Wait(request, status);
505 _doEventLeave(eMpiName::Wait);
514waitAll(
int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[])
516 _doEventEnter(eMpiName::Waitall);
517 int r = MPI_Waitall(count, array_of_requests, array_of_statuses);
518 _doEventLeave(eMpiName::Waitall);
527testSome(
int incount, MPI_Request array_of_requests[],
int *outcount,
528 int array_of_indices[], MPI_Status array_of_statuses[])
530 _doEventEnter(eMpiName::Testsome);
531 int r = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
532 _doEventLeave(eMpiName::Testsome);
541waitSome(
int incount, MPI_Request array_of_requests[],
int *outcount,
542 int array_of_indices[], MPI_Status array_of_statuses[])
544 _doEventEnter(eMpiName::Waitsome);
545 int r = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
546 _doEventLeave(eMpiName::Waitsome);
553void Otf2MpiProfiling::
554_doEventEnter(eMpiName event_name)
557 static_cast<uint32_t
>(event_name));
563void Otf2MpiProfiling::
564_doEventLeave(eMpiName event_name)
567 static_cast<uint32_t
>(event_name));
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 -*-