15#include "arccore/message_passing_mpi/internal/MessagePassingMpiEnum.h"
16#include "arcane/utils/FatalErrorException.h"
17#include "arcane/std/internal/Otf2MpiProfiling.h"
24using namespace MessagePassing::Mpi;
32 using ReturnType = Otf2MpiProfiling::ReturnType;
35 getSizeOfMpiType(MPI_Datatype datatype)
38 MPI_Type_size(datatype, &tmp_size);
39 return static_cast<uint64_t
>(tmp_size);
53: m_otf2_wrapper(otf2_wrapper)
62broadcast(
void* buffer,
int count, MPI_Datatype datatype,
int root, MPI_Comm comm)
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);
68 _doEventEnter(eMpiName::Bcast);
72 int r = MPI_Bcast(buffer, count, datatype, root, comm);
75 OTF2_COLLECTIVE_OP_BCAST, 0 , root,
76 byte_send , byte_recv );
78 _doEventLeave(eMpiName::Bcast);
87gather(
const void* sendbuf,
int sendcount, MPI_Datatype sendtype,
void* recvbuf,
88 int recvcount, MPI_Datatype recvtype,
int root, MPI_Comm comm)
91 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
92 const uint64_t byte_recv(m_otf2_wrapper->getMpiRank() == root ? getSizeOfMpiType(recvtype) * recvcount : 0);
94 _doEventEnter(eMpiName::Gather);
98 int r = MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
101 OTF2_COLLECTIVE_OP_GATHER, 0 , root,
102 byte_send , byte_recv );
104 _doEventLeave(eMpiName::Gather);
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)
117 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
118 uint64_t byte_recv(0);
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)
125 _doEventEnter(eMpiName::Gatherv);
129 int r = MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
132 OTF2_COLLECTIVE_OP_GATHERV, 0 , root,
133 byte_send , byte_recv );
135 _doEventLeave(eMpiName::Gatherv);
144allGather(
const void* sendbuf,
int sendcount, MPI_Datatype sendtype,
void* recvbuf,
145 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
148 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
149 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount * m_otf2_wrapper->getMpiNbRank());
151 _doEventEnter(eMpiName::Allgather);
155 int r = MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
158 OTF2_COLLECTIVE_OP_ALLGATHER, 0 , OTF2_UNDEFINED_UINT32 ,
159 byte_send , byte_recv );
161 _doEventLeave(eMpiName::Allgather);
170allGatherVariable(
const void* sendbuf,
int sendcount, MPI_Datatype sendtype,
void* recvbuf,
171 const int* recvcounts,
const int* displs, MPI_Datatype recvtype, MPI_Comm comm)
174 const uint64_t byte_send(getSizeOfMpiType(sendtype) * sendcount);
175 uint64_t byte_recv(0);
177 for (
int i(0); i < m_otf2_wrapper->getMpiNbRank(); ++i)
178 byte_recv += recvcounts[i] * getSizeOfMpiType(recvtype);
180 _doEventEnter(eMpiName::Allgatherv);
184 int r = MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
187 OTF2_COLLECTIVE_OP_ALLGATHERV, 0 , OTF2_UNDEFINED_UINT32 ,
188 byte_send , byte_recv );
190 _doEventLeave(eMpiName::Allgatherv);
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)
204 const uint64_t byte_recv(getSizeOfMpiType(recvtype) * recvcount);
205 uint64_t byte_send(0);
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)
213 _doEventEnter(eMpiName::Scatterv);
217 int r = MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
220 OTF2_COLLECTIVE_OP_SCATTERV, 0 , root,
221 byte_send , byte_recv );
223 _doEventLeave(eMpiName::Scatterv);
232allToAll(
const void* sendbuf,
int sendcount, MPI_Datatype sendtype,
void* recvbuf,
233 int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
235 _doEventEnter(eMpiName::Alltoall);
239 int r = MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
242 OTF2_COLLECTIVE_OP_ALLTOALL, 0 , OTF2_UNDEFINED_UINT32 ,
245 _doEventLeave(eMpiName::Alltoall);
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)
258 _doEventEnter(eMpiName::Alltoallv);
262 int r = MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf, recvcounts, rdispls, recvtype, comm);
265 OTF2_COLLECTIVE_OP_ALLTOALLV, 0 , OTF2_UNDEFINED_UINT32 ,
268 _doEventLeave(eMpiName::Alltoallv);
279 _doEventEnter(eMpiName::Barrier);
283 int r = MPI_Barrier(comm);
286 OTF2_COLLECTIVE_OP_BARRIER, 0 , OTF2_UNDEFINED_UINT32 ,
289 _doEventLeave(eMpiName::Barrier);
298reduce(
const void* sendbuf,
void* recvbuf,
int count, MPI_Datatype datatype,
299 MPI_Op op,
int root, MPI_Comm comm)
301 _doEventEnter(eMpiName::Reduce);
305 int r = MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
308 OTF2_COLLECTIVE_OP_REDUCE, 0 , root,
311 _doEventLeave(eMpiName::Reduce);
320allReduce(
const void* sendbuf,
void* recvbuf,
int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
323 const uint64_t byte_send(getSizeOfMpiType(datatype) * count);
324 const uint64_t byte_recv(getSizeOfMpiType(datatype) * count);
326 _doEventEnter(eMpiName::Allreduce);
330 int r = MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
333 OTF2_COLLECTIVE_OP_ALLREDUCE, 0 , OTF2_UNDEFINED_UINT32 ,
334 byte_send , byte_recv );
336 _doEventLeave(eMpiName::Allreduce);
345scan(
const void* sendbuf,
void* recvbuf,
int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
347 _doEventEnter(eMpiName::Scan);
351 int r = MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
354 MPI_Type_size(datatype, &type_size);
355 const uint64_t size =
static_cast<uint64_t
>(type_size) *
static_cast<uint64_t
>(count);
357 OTF2_COLLECTIVE_OP_SCAN, 0 , OTF2_UNDEFINED_UINT32 ,
360 _doEventLeave(eMpiName::Scan);
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)
373 _doEventEnter(eMpiName::Sendrecv);
375 int r = MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype,
376 source, recvtag, comm, status);
378 _doEventLeave(eMpiName::Sendrecv);
387iSend(
const void* buf,
int count, MPI_Datatype datatype,
int dest,
int tag,
388 MPI_Comm comm, MPI_Request* request)
390 _doEventEnter(eMpiName::Isend);
392 int r = MPI_Isend(buf, count, datatype, dest, tag, comm, request);
395 uint64_t request_id = MPI_Request_c2f(*request);
397 dest, 0 , tag, getSizeOfMpiType(datatype) * count, request_id);
399 _doEventLeave(eMpiName::Isend);
408send(
const void* buf,
int count, MPI_Datatype datatype,
int dest,
int tag, MPI_Comm comm)
410 _doEventEnter(eMpiName::Send);
412 int r = MPI_Send(buf, count, datatype, dest, tag, comm);
415 dest, 0 , tag, getSizeOfMpiType(datatype) * count);
417 _doEventLeave(eMpiName::Send);
426iRecv(
void* buf,
int count, MPI_Datatype datatype,
int source,
int tag,
427 MPI_Comm comm, MPI_Request* request)
429 _doEventEnter(eMpiName::Irecv);
431 int r = MPI_Irecv(buf, count, datatype, source, tag, comm, request);
436 _doEventLeave(eMpiName::Irecv);
445recv(
void* buf,
int count, MPI_Datatype datatype,
int source,
int tag, MPI_Comm comm, MPI_Status* status)
447 _doEventEnter(eMpiName::Recv);
448 int r = MPI_Recv(buf, count, datatype, source, tag, comm, status);
451 source, 0 , tag, getSizeOfMpiType(datatype) * count);
453 _doEventLeave(eMpiName::Recv);
462test(MPI_Request* request,
int* flag, MPI_Status* status)
464 _doEventEnter(eMpiName::Test);
465 if (!MPI_Test(request, flag, status)) {
466 uint64_t request_id = MPI_Request_c2f(*request);
469 _doEventLeave(eMpiName::Test);
478probe(
int source,
int tag, MPI_Comm comm, MPI_Status* status)
480 _doEventEnter(eMpiName::Probe);
481 int r = MPI_Probe(source, tag, comm, status);
482 _doEventLeave(eMpiName::Probe);
491getCount(
const MPI_Status* status, MPI_Datatype datatype,
int* count)
493 _doEventEnter(eMpiName::Get_count);
494 int r = MPI_Get_count(status, datatype, count);
495 _doEventLeave(eMpiName::Get_count);
504wait(MPI_Request* request, MPI_Status* status)
506 _doEventEnter(eMpiName::Wait);
507 int r = MPI_Wait(request, status);
508 _doEventLeave(eMpiName::Wait);
517waitAll(
int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[])
519 _doEventEnter(eMpiName::Waitall);
520 int r = MPI_Waitall(count, array_of_requests, array_of_statuses);
521 _doEventLeave(eMpiName::Waitall);
530testSome(
int incount, MPI_Request array_of_requests[],
int* outcount,
531 int array_of_indices[], MPI_Status array_of_statuses[])
533 _doEventEnter(eMpiName::Testsome);
534 int r = MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
535 _doEventLeave(eMpiName::Testsome);
544waitSome(
int incount, MPI_Request array_of_requests[],
int* outcount,
545 int array_of_indices[], MPI_Status array_of_statuses[])
547 _doEventEnter(eMpiName::Waitsome);
548 int r = MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_statuses);
549 _doEventLeave(eMpiName::Waitsome);
556void Otf2MpiProfiling::
557_doEventEnter(eMpiName event_name)
560 static_cast<uint32_t
>(event_name));
566void Otf2MpiProfiling::
567_doEventLeave(eMpiName event_name)
570 static_cast<uint32_t
>(event_name));
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 --