Salome HOME
ParaUMesh.redistributeCells implementation.
[tools/medcoupling.git] / src / ParaMEDMEM / CommInterface.hxx
1 // Copyright (C) 2007-2020  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 #pragma once
21
22 #include "ParaIdType.hxx"
23 #include "MEDCouplingMemArray.hxx"
24
25 #include <mpi.h>
26
27 #include <memory>
28 #include <numeric>
29
30 namespace MEDCoupling
31 {
32   template<class T>
33   struct ParaTraits
34   {
35     using EltType = T;
36   };
37   
38   template<>
39   struct ParaTraits<double>
40   {
41     static MPI_Datatype MPIDataType;
42   };
43
44   template<>
45   struct ParaTraits<Int32>
46   {
47     static MPI_Datatype MPIDataType;
48   };
49
50   template<>
51   struct ParaTraits<Int64>
52   {
53     static MPI_Datatype MPIDataType;
54   };
55
56   class CommInterface
57   {
58   public:
59     CommInterface() { }
60     virtual ~CommInterface() { }
61     int worldSize() const {
62       int size;
63       MPI_Comm_size(MPI_COMM_WORLD, &size);
64       return size;}
65     int commSize(MPI_Comm comm, int* size) const { return MPI_Comm_size(comm,size); }
66     int commRank(MPI_Comm comm, int* rank) const { return MPI_Comm_rank(comm,rank); }
67     int commGroup(MPI_Comm comm, MPI_Group* group) const { return MPI_Comm_group(comm, group); }
68     int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const { return MPI_Group_incl(group, size, ranks, group_output); }
69     int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const { return MPI_Comm_create(comm,group,comm_output); }
70     int groupFree(MPI_Group* group) const { return MPI_Group_free(group); }
71     int commFree(MPI_Comm* comm) const { return MPI_Comm_free(comm); }
72
73     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); }
74     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); }
75     int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
76                  int dest, int sendtag, void* recvbuf, int recvcount, 
77                  MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
78                  MPI_Status* status) { return MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype, source, recvtag, comm,status); }
79
80     int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
81               int tag, MPI_Comm comm, MPI_Request *request) const { return MPI_Isend(buffer,count, datatype, target, tag, comm, request); }
82     int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
83               int tag, MPI_Comm comm, MPI_Request* request) const { return MPI_Irecv(buffer,count, datatype, source, tag, comm, request); }
84
85     int wait(MPI_Request *request, MPI_Status *status) const { return MPI_Wait(request, status); }
86     int test(MPI_Request *request, int *flag, MPI_Status *status) const { return MPI_Test(request, flag, status); }
87     int requestFree(MPI_Request *request) const { return MPI_Request_free(request); }
88     int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const { return MPI_Waitany(count, array_of_requests, index, status); }
89     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); }
90     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); }
91     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); }
92     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); }
93     int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
94                  int *array_of_indices, MPI_Status *array_of_status) const { return MPI_Testsome(incount, array_of_requests, outcount, array_of_indices, array_of_status); }
95     int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const { return MPI_Probe(source, tag, comm, status) ; }
96     int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const { return MPI_Iprobe(source, tag, comm, flag, status) ; }
97     int cancel(MPI_Request *request) const { return MPI_Cancel(request); }
98     int testCancelled(MPI_Status *status, int *flag) const { return MPI_Test_cancelled(status, flag); }
99     int barrier(MPI_Comm comm) const { return MPI_Barrier(comm); }
100     int errorString(int errorcode, char *string, int *resultlen) const { return MPI_Error_string(errorcode, string, resultlen); }
101     int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const { return MPI_Get_count(status, datatype, count); }
102
103     int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const { return MPI_Bcast(buffer, count,  datatype, root, comm); }
104     int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
105                   void* recvbuf, int recvcount, MPI_Datatype recvtype,
106                   MPI_Comm comm) const { return MPI_Allgather(sendbuf,sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
107     int allGatherV(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[],
108                    const int displs[], MPI_Datatype recvtype, MPI_Comm comm) const { return MPI_Allgatherv(sendbuf,sendcount,sendtype,recvbuf,recvcounts,displs,recvtype,comm); }
109     int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
110                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
111                  MPI_Comm comm) const { return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
112     int allToAllV(const void* sendbuf, int* sendcounts, int* senddispls,
113                   MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
114                   int* recvdispls, MPI_Datatype recvtype, 
115                   MPI_Comm comm) const { return MPI_Alltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm); }
116
117     int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype,
118                MPI_Op op, int root, MPI_Comm comm) const { return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm); }
119     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); }
120   public:
121     void allGatherArrays(MPI_Comm comm, const DataArrayIdType *array, std::vector< MCAuto<DataArrayIdType> >& arraysOut) const;
122     int allGatherArrays(MPI_Comm comm, const DataArrayIdType *array, std::unique_ptr<mcIdType[]>& result, std::unique_ptr<mcIdType[]>& resultIndex) const;
123     void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayIdType> >& arrays, std::vector< MCAuto<DataArrayIdType> >& arraysOut) const;
124     void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayDouble> >& arrays, std::vector< MCAuto<DataArrayDouble> >& arraysOut) const;
125     void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayDouble> >& arrays, MCAuto<DataArrayDouble>& arraysOut) const;
126     template<class T>
127     int allToAllArraysT2(MPI_Comm comm, const std::vector< MCAuto<typename Traits<T>::ArrayType> >& arrays, MCAuto<typename Traits<T>::ArrayType>& arrayOut, std::unique_ptr<mcIdType[]>& nbOfElems2, mcIdType& nbOfComponents) const
128     {
129       using DataArrayT = typename Traits<T>::ArrayType;
130       int size;
131       this->commSize(comm,&size);
132       if( arrays.size() != ToSizeT(size) )
133         throw INTERP_KERNEL::Exception("AllToAllArrays : internal error ! Invalid size of input array.");
134         
135       std::vector< const DataArrayT *> arraysBis(FromVecAutoToVecOfConst<DataArrayT>(arrays));
136       std::unique_ptr<mcIdType[]> nbOfElems3(new mcIdType[size]);
137       nbOfElems2.reset(new mcIdType[size]);
138       nbOfComponents = std::numeric_limits<mcIdType>::max();
139       for(int curRk = 0 ; curRk < size ; ++curRk)
140       {
141         mcIdType curNbOfCompo( ToIdType( arrays[curRk]->getNumberOfComponents() ) );
142         if(nbOfComponents != std::numeric_limits<mcIdType>::max())
143         {
144           if( nbOfComponents != curNbOfCompo )
145             throw INTERP_KERNEL::Exception("AllToAllArrays : internal error ! Nb of components is not homogeneous !");
146         }
147         else
148         {
149           nbOfComponents = curNbOfCompo;
150         }
151         nbOfElems3[curRk] = arrays[curRk]->getNbOfElems();
152       }
153       this->allToAll(nbOfElems3.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
154       mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0));
155       arrayOut = DataArrayT::New();
156       arrayOut->alloc(nbOfCellIdsSum/nbOfComponents,nbOfComponents);
157       std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems3,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
158       std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) );
159       {
160         MCAuto<DataArrayT> arraysAcc(DataArrayT::Aggregate(arraysBis));
161         this->allToAllV(arraysAcc->begin(),nbOfElemsInt.get(),offsetsIn.get(),ParaTraits<T>::MPIDataType,
162                         arrayOut->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),ParaTraits<T>::MPIDataType,comm);
163       }
164       return size;
165     }
166
167     template<class T>
168     void allToAllArraysT(MPI_Comm comm, const std::vector< MCAuto<typename Traits<T>::ArrayType> >& arrays, std::vector< MCAuto<typename Traits<T>::ArrayType> >& arraysOut) const
169     {
170       using DataArrayT = typename Traits<T>::ArrayType;
171       MCAuto<DataArrayT> cellIdsFromProcs;
172       std::unique_ptr<mcIdType[]> nbOfElems2;
173       mcIdType nbOfComponents(0);
174       int size(this->allToAllArraysT2<T>(comm,arrays,cellIdsFromProcs,nbOfElems2,nbOfComponents));
175       std::unique_ptr<mcIdType[]> offsetsOutIdType( CommInterface::ComputeOffset(nbOfElems2,size) );
176       // build output arraysOut by spliting cellIdsFromProcs into parts
177       arraysOut.resize(size);
178       for(int curRk = 0 ; curRk < size ; ++curRk)
179       {
180         arraysOut[curRk] = DataArrayT::NewFromArray(cellIdsFromProcs->begin()+offsetsOutIdType[curRk],cellIdsFromProcs->begin()+offsetsOutIdType[curRk]+nbOfElems2[curRk]);
181         arraysOut[curRk]->rearrange(nbOfComponents);
182       }
183     }
184   public:
185
186     /*!
187     * \a counts is expected to be an array of array length. This method returns an array of split array.
188     */
189     static std::unique_ptr<mcIdType[]> SplitArrayOfLength(const std::unique_ptr<mcIdType[]>& counts, std::size_t countsSz, int rk, int size)
190     {
191       std::unique_ptr<mcIdType[]> ret(new mcIdType[countsSz]);
192       for(std::size_t i=0;i<countsSz;++i)
193       {
194         mcIdType a,b;
195         DataArray::GetSlice(0,counts[i],1,rk,size,a,b);
196         ret[i] = b-a;
197       }
198       return ret;
199     }
200
201     /*!
202     * Helper of alltoallv and allgatherv
203     */
204     template<class T>
205     static std::unique_ptr<int []> ToIntArray(const std::unique_ptr<T []>& arr, std::size_t size)
206     {
207       std::unique_ptr<int []> ret(new int[size]);
208       std::copy(arr.get(),arr.get()+size,ret.get());
209       return ret;
210     }
211     
212     /*!
213     * Helper of alltoallv and allgatherv
214     */
215     template<class T>
216     static std::unique_ptr<T []> ComputeOffset(const std::unique_ptr<T []>& counts, std::size_t sizeOfCounts)
217     {
218       std::unique_ptr<T []> ret(new T[sizeOfCounts]);
219       ret[0] = static_cast<T>(0);
220       for(std::size_t i = 1 ; i < sizeOfCounts ; ++i)
221       {
222         ret[i] = ret[i-1] + counts[i-1];
223       }
224       return ret;
225     }
226
227     /*!
228     * Helper of alltoallv and allgatherv
229     */
230     template<class T>
231     static std::unique_ptr<T []> ComputeOffsetFull(const std::unique_ptr<T []>& counts, std::size_t sizeOfCounts)
232     {
233       std::unique_ptr<T []> ret(new T[sizeOfCounts+1]);
234       ret[0] = static_cast<T>(0);
235       for(std::size_t i = 1 ; i < sizeOfCounts+1 ; ++i)
236       {
237         ret[i] = ret[i-1] + counts[i-1];
238       }
239       return ret;
240     }
241   };
242 }