]> SALOME platform Git repositories - tools/medcoupling.git/blob - src/ParaMEDMEM/CommInterface.hxx
Salome HOME
Indices are stored as mcIdType type instead of int to support switch to 64bits indexing
[tools/medcoupling.git] / src / ParaMEDMEM / CommInterface.hxx
1 // Copyright (C) 2007-2019  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #ifndef __COMMINTERFACE_HXX__
21 #define __COMMINTERFACE_HXX__
22
23 #include <mpi.h>
24
25 #include "ParaIdType.hxx"
26
27 namespace MEDCoupling
28 {
29
30   class CommInterface
31   {
32   public:
33     CommInterface(){}
34     virtual ~CommInterface(){}
35     int worldSize() const {
36       int size;
37       MPI_Comm_size(MPI_COMM_WORLD, &size);
38       return size;}
39     int commSize(MPI_Comm comm, int* size) const { return MPI_Comm_size(comm,size); }
40     int commRank(MPI_Comm comm, int* rank) const { return MPI_Comm_rank(comm,rank); }
41     int commGroup(MPI_Comm comm, MPI_Group* group) const { return MPI_Comm_group(comm, group); }
42     int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const { return MPI_Group_incl(group, size, ranks, group_output); }
43     int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const { return MPI_Comm_create(comm,group,comm_output); }
44     int groupFree(MPI_Group* group) const { return MPI_Group_free(group); }
45     int commFree(MPI_Comm* comm) const { return MPI_Comm_free(comm); }
46
47     int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const { return MPI_Send(buffer,count, datatype, target, tag, comm); }
48     int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const { return MPI_Recv(buffer,count, datatype, source, tag, comm, status); }
49     int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
50                  int dest, int sendtag, void* recvbuf, int recvcount, 
51                  MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
52                  MPI_Status* status) { return MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype, source, recvtag, comm,status); }
53
54     int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
55               int tag, MPI_Comm comm, MPI_Request *request) const { return MPI_Isend(buffer,count, datatype, target, tag, comm, request); }
56     int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
57               int tag, MPI_Comm comm, MPI_Request* request) const { return MPI_Irecv(buffer,count, datatype, source, tag, comm, request); }
58
59     int wait(MPI_Request *request, MPI_Status *status) const { return MPI_Wait(request, status); }
60     int test(MPI_Request *request, int *flag, MPI_Status *status) const { return MPI_Test(request, flag, status); }
61     int requestFree(MPI_Request *request) const { return MPI_Request_free(request); }
62     int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const { return MPI_Waitany(count, array_of_requests, index, status); }
63     int testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status) const { return MPI_Testany(count, array_of_requests, index, flag, status); }
64     int waitall(int count, MPI_Request *array_of_requests, MPI_Status *array_of_status) const { return MPI_Waitall(count, array_of_requests, array_of_status); }
65     int testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_status) const { return MPI_Testall(count, array_of_requests, flag, array_of_status); }
66     int waitsome(int incount, MPI_Request *array_of_requests,int *outcount, int *array_of_indices, MPI_Status *array_of_status) const { return MPI_Waitsome(incount, array_of_requests, outcount, array_of_indices, array_of_status); }
67     int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
68                  int *array_of_indices, MPI_Status *array_of_status) const { return MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_status); }
69     int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const { return MPI_Probe(source, tag, comm, status) ; }
70     int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const { return MPI_Iprobe(source, tag, comm, flag, status) ; }
71     int cancel(MPI_Request *request) const { return MPI_Cancel(request); }
72     int testCancelled(MPI_Status *status, int *flag) const { return MPI_Test_cancelled(status, flag); }
73     int barrier(MPI_Comm comm) const { return MPI_Barrier(comm); }
74     int errorString(int errorcode, char *string, int *resultlen) const { return MPI_Error_string(errorcode, string, resultlen); }
75     int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const { return MPI_Get_count(status, datatype, count); }
76
77     int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const { return MPI_Bcast(buffer, count,  datatype, root, comm); }
78     int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
79                   void* recvbuf, int recvcount, MPI_Datatype recvtype,
80                   MPI_Comm comm) const { return MPI_Allgather(sendbuf,sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
81     int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
82                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
83                  MPI_Comm comm) const { return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
84     int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
85                   MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
86                   int* recvdispls, MPI_Datatype recvtype, 
87                   MPI_Comm comm) const { return MPI_Alltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm); }
88
89     int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype,
90                MPI_Op op, int root, MPI_Comm comm) const { return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm); }
91     int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const { return MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm); }
92   };
93 }
94
95 #endif /*COMMINTERFACE_HXX_*/