source: XIOS/dev/branch_yushan/extern/remap/src/mpi_routing.cpp @ 1053

Last change on this file since 1053 was 1053, checked in by yushan, 7 years ago

ep_lib namespace specified when netcdf involved

File size: 18.3 KB
Line 
1#include "mpi_routing.hpp"
2#include "mpi.hpp"
3#include "node.hpp"
4#include "elt.hpp"
5#include "timerRemap.hpp"
6#include <iostream>
7#ifdef _usingEP
8#include "ep_declaration.hpp"
9#endif
10
11namespace sphereRemap {
12
13const int verbose = 0;
14
15CMPIRouting::CMPIRouting(MPI_Comm comm) : communicator(comm)
16{
17        MPI_Comm_rank(comm, &mpiRank);
18        MPI_Comm_size(comm, &mpiSize);
19}
20
21/* sparse alltoallv when it is known that only a subset of ranks will communicate,
22    but message lengths are *known* to receiver */
23template <typename T>
24void alltoalls_known(const vector<vector<T> >& send, vector<vector<T> >& recv, const vector<int>& ranks, MPI_Comm communicator)
25{
26        vector<MPI_Request> request(ranks.size() * 2);
27        vector<MPI_Status>  status(ranks.size() * 2);
28
29        // communicate data
30        int nbRequest = 0;
31        for (int i = 0; i < ranks.size(); i++)
32                if (recv[i].size())
33                        MPI_Irecv(&recv[i][0], recv[i].size()*sizeof(T), MPI_CHAR, ranks[i], 0, communicator, &request[nbRequest++]);
34        for (int i = 0; i < ranks.size(); i++)
35                if (send[i].size())
36                        MPI_Isend((void *) &send[i][0], send[i].size()*sizeof(T), MPI_CHAR, ranks[i], 0, communicator, &request[nbRequest++]);
37        MPI_Waitall(nbRequest, &request[0], &status[0]);
38}
39
40/* sparse alltoallv when it is known that only a subset of ranks will communicate,
41    but message lengths are *unknown* to receiver */
42template <typename T>
43void alltoalls_unknown(const vector<vector<T> >& send, vector<vector<T> >& recv, const vector<int>& ranks, MPI_Comm communicator)
44{
45        vector<MPI_Request> request(ranks.size() * 2);
46        vector<MPI_Status>  status(ranks.size() * 2);
47
48        // communicate sizes
49        int nbRequest = 0;
50        vector<int> sendSizes(ranks.size());
51        vector<int> recvSizes(ranks.size());
52        for (int i = 0; i < ranks.size(); i++)
53                sendSizes[i] = send[i].size();
54        for (int i = 0; i < ranks.size(); i++)
55                MPI_Irecv(&recvSizes[i], 1, MPI_INT, ranks[i], 0, communicator, &request[nbRequest++]);
56        for (int i = 0; i < ranks.size(); i++)
57                MPI_Isend(&sendSizes[i], 1, MPI_INT, ranks[i], 0, communicator, &request[nbRequest++]);
58        MPI_Waitall(nbRequest, &request[0], &status[0]);
59
60        // allocate
61        for (int i = 0; i < ranks.size(); i++)
62                if (recvSizes[i]) recv[i].resize(recvSizes[i]);
63        // communicate data
64        alltoalls_known(send, recv, ranks, communicator);
65}
66
67void setElementsSendCnt(int route, vector<int>& sendCnt)
68{
69        sendCnt[route]++;
70}
71
72void setElementsSendCnt(const vector<int>& route, vector<int>& sendCnt)
73{
74        for (int i = 0; i < route.size(); i++)
75                sendCnt[route[i]]++;
76}
77
78void setTargetElementIndex(int route, vector<int>& indices, const int *rankToIndex)
79{
80        int index = rankToIndex[route];
81        indices.push_back(index);
82}
83
84void setTargetElementIndex(const vector<int>& route, vector<int>& indices, const int *rankToIndex)
85{
86        for (int i = 0; i < route.size(); i++)
87        {
88                int index = rankToIndex[route[i]];
89                indices.push_back(index);
90        }
91}
92
93template<typename T>
94void CMPIRouting::init(const vector<T>& route, CMPICascade *cascade)
95{
96        vector<int> nbElementSend(mpiSize);
97        int *toSend = new int[mpiSize];
98        int *recvCount = new int[mpiSize];
99        int *targetRankToIndex;
100
101        for (int i = 0; i < route.size(); i++)
102                setElementsSendCnt(route[i], nbElementSend);
103
104        nbTarget = 0;
105vector<int> destRanks;
106        for (int i = 0; i < mpiSize; i++)
107        {
108                if (nbElementSend[i] == 0)
109                        toSend[i] = 0;
110                else
111                {
112destRanks.push_back(i);
113                        toSend[i] = 1;
114                        nbTarget++;
115                }
116                recvCount[i] = 1;
117        }
118//int recvCntDeb = (stree) ? stree->scatter_reduce(destRanks) : -1;
119
120        MPI_Barrier(communicator);
121        CTimer::get("CMPIRouting::init(reduce_scatter)").reset();
122        CTimer::get("CMPIRouting::init(reduce_scatter)").resume();
123        MPI_Reduce_scatter(toSend, &nbSource, recvCount, MPI_INT, MPI_SUM, communicator);
124        CTimer::get("CMPIRouting::init(reduce_scatter)").suspend();
125        CTimer::get("CMPIRouting::init(reduce_scatter)").print();
126
127        MPI_Info info_null;
128
129        MPI_Alloc_mem(nbTarget *sizeof(int), info_null, &targetRank);
130        MPI_Alloc_mem(nbSource *sizeof(int), info_null, &sourceRank);
131
132        targetRankToIndex = new int[mpiSize];
133        int index = 0;
134        for (int i = 0; i < mpiSize; i++)
135        {
136                if (toSend[i] == 1)
137                {
138                        targetRank[index] = i;
139                        targetRankToIndex[i] = index;
140                        index++;
141                }
142        }
143
144        MPI_Barrier(communicator);
145        CTimer::get("CMPIRouting::init(get_source)").reset();
146        CTimer::get("CMPIRouting::init(get_source)").resume();
147
148        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
149        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
150
151        int indexRequest = 0;
152        if (verbose) cout << "CMPIRouting::init     nbSource : " << nbSource << " nbTarget : " << nbTarget << endl;
153//      assert(recvCntDeb == -1 || recvCntDeb == nbSource);
154//cout << mpiRank <<  "DEB : " << recvCntDeb << "    nbSource " << nbSource << " nbTarget : " << nbTarget << endl;
155        for (int i = 0; i < nbSource; i++)
156        {
157                MPI_Irecv(&sourceRank[i], 1, MPI_INT, MPI_ANY_SOURCE, 0, communicator, &request[indexRequest]);
158                indexRequest++;
159        }
160        MPI_Barrier(communicator);
161        for (int i = 0; i < nbTarget; i++)
162        {
163                MPI_Isend(&mpiRank, 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
164                indexRequest++;
165        }
166        MPI_Waitall(indexRequest, request, status);
167        MPI_Barrier(communicator);  //needed
168        CTimer::get("CMPIRouting::init(get_source)").suspend();
169        CTimer::get("CMPIRouting::init(get_source)").print();
170
171        CTimer::get("CMPIRouting::init(get_source)").reset();
172        CTimer::get("CMPIRouting::init(get_source)").resume();
173
174        indexRequest = 0;
175        for (int i = 0; i < nbSource; i++)
176        {
177                MPI_Irecv(&sourceRank[i], 1, MPI_INT, MPI_ANY_SOURCE, 0, communicator, &request[indexRequest]);
178                indexRequest++;
179        }
180
181        for (int i = 0; i < nbTarget; i++)
182        {
183                MPI_Isend(&mpiRank, 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
184                indexRequest++;
185        }
186        MPI_Waitall(indexRequest, request, status);
187        MPI_Barrier(communicator);
188        CTimer::get("CMPIRouting::init(get_source)").suspend();
189        CTimer::get("CMPIRouting::init(get_source)").print();
190
191        CTimer::get("CMPIRouting::init(send_element)").reset();
192        CTimer::get("CMPIRouting::init(send_element)").resume();
193
194        nbTargetElement.resize(nbTarget);
195        nbSourceElement.resize(nbSource);
196
197        for (int i = 0; i < route.size(); i++)
198                setTargetElementIndex(route[i], targetElementIndex, targetRankToIndex);
199
200        for (int i = 0; i < targetElementIndex.size(); i++)
201                nbTargetElement[targetElementIndex[i]]++;
202
203        indexRequest = 0;
204        totalSourceElement = 0;
205        totalTargetElement = 0;
206        for (int i = 0; i < nbSource; i++)
207        {
208                MPI_Irecv(&nbSourceElement[i], 1, MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
209                indexRequest++;
210        }
211
212        for (int i = 0; i < nbTarget; i++)
213        {
214                totalTargetElement += nbTargetElement[i];
215                MPI_Isend(&nbTargetElement[i], 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
216                indexRequest++;
217        }
218
219        MPI_Waitall(indexRequest, request, status);
220
221        CTimer::get("CMPIRouting::init(send_element)").suspend();
222        CTimer::get("CMPIRouting::init(send_element)").print();
223
224        totalSourceElement = 0;
225        for (int i = 0; i < nbSource; i++)
226                totalSourceElement += nbSourceElement[i];
227
228        sourceElementIndex.resize(totalSourceElement);
229
230        totalSourceElement = 0;
231        for (int i = 0; i < nbSource; i++)
232        {
233                for (int j = 0; j < nbSourceElement[i]; j++)
234                {
235                        sourceElementIndex[totalSourceElement] = i;
236                        totalSourceElement++;
237                }
238        }
239
240        delete[]  toSend;
241        delete[]  recvCount;
242        delete[]  request;
243        delete[]  status;
244}
245
246
247int CMPIRouting::getTotalSourceElement(void)
248{
249        return totalSourceElement;
250}
251
252template<typename T>
253void CMPIRouting::transferToTarget(T* targetElements, T* sourceElements)
254{
255        char** targetBuffer = new char*[nbTarget];
256        int* indexTargetBuffer = new int[nbTarget];
257
258        for(int i=0;i<nbTarget; i++)
259        {
260                targetBuffer[i] = new char[sizeof(T)*nbTargetElement[i]];
261                indexTargetBuffer[i]= 0;
262        }
263
264        char** sourceBuffer = new char*[nbSource];
265        int* indexSourceBuffer = new int[nbSource];
266
267        for(int i=0;i<nbSource; i++)
268        {
269                sourceBuffer[i] = new char[sizeof(T)*nbSourceElement[i]];
270                indexSourceBuffer[i]= 0;
271        }
272
273        // pack the data
274        int index;
275        for(int i=0;i<totalTargetElement;i++)
276        {
277                index=targetElementIndex[i];
278                *((T*) &(targetBuffer[index][indexTargetBuffer[index]])) = targetElements[i];
279                indexTargetBuffer[index]+=sizeof(T);
280        }
281
282
283        MPI_Request* request=new MPI_Request[nbSource+nbTarget];
284        MPI_Status*  status=new MPI_Status[nbSource+nbTarget];
285        int indexRequest=0;
286
287        MPI_Barrier(communicator);
288        CTimer::get("CMPIRouting::transferToTarget").reset();
289        CTimer::get("CMPIRouting::transferToTarget").resume();
290
291        for(int i=0; i<nbSource; i++)
292        {
293                MPI_Irecv(sourceBuffer[i],nbSourceElement[i]*sizeof(T),MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
294                indexRequest++;
295        }
296
297        for(int i=0;i<nbTarget; i++)
298        {
299                MPI_Isend(targetBuffer[i],nbTargetElement[i]*sizeof(T), MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
300                indexRequest++;
301        }
302
303        MPI_Waitall(indexRequest,request,status);
304
305        CTimer::get("CMPIRouting::transferToTarget").suspend();
306        CTimer::get("CMPIRouting::transferToTarget").print();
307        MPI_Barrier(communicator);
308
309        // unpack the data
310        for(int i=0;i<totalSourceElement;i++)
311        {
312                index=sourceElementIndex[i];
313                sourceElements[i] = *((T*) &(sourceBuffer[index][indexSourceBuffer[index]]));
314                indexSourceBuffer[index]+=sizeof(T);
315        }
316
317        for(int i=0;i<nbTarget; i++) delete[] targetBuffer[i];
318        for(int i=0;i<nbSource; i++) delete[] sourceBuffer[i];
319        delete[] targetBuffer;
320        delete[] indexTargetBuffer;
321        delete[] sourceBuffer;
322        delete[] indexSourceBuffer;
323        delete[] request;
324        delete[] status;
325}
326
327
328template<typename T, typename t_pack, typename t_unpack>
329void CMPIRouting::transferToTarget(T* targetElements, T* sourceElements, t_pack pack, t_unpack unpack)
330{
331        char** targetBuffer = new char*[nbTarget];
332        int* indexTargetBuffer = new int[nbTarget];
333        int* targetMessageSize = new int[nbTarget];
334        int* sourceMessageSize = new int[nbSource];
335        int index;
336
337        // compute target buffer size
338        for (int i = 0; i < nbTarget; i++)
339                targetMessageSize[i] = 0;
340
341        for (int i = 0; i < totalTargetElement; i++)
342        {
343                index = targetElementIndex[i];
344                pack(targetElements[i], NULL, targetMessageSize[index]);
345        }
346
347        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
348        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
349        int indexRequest = 0;
350
351        MPI_Barrier(communicator);
352        CTimer::get("CMPIRouting::transferToTarget(messageSize)").reset();
353        CTimer::get("CMPIRouting::transferToTarget(messageSize)").resume();
354
355        for(int i=0; i<nbSource; i++)
356        {
357                MPI_Irecv(&sourceMessageSize[i],1,MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
358                indexRequest++;
359        }
360
361        for(int i=0; i<nbTarget; i++)
362        {
363                MPI_Isend(&targetMessageSize[i],1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
364                indexRequest++;
365        }
366
367        MPI_Waitall(indexRequest,request,status);
368
369        MPI_Barrier(communicator);
370        CTimer::get("CMPIRouting::transferToTarget(messageSize)").suspend();
371        CTimer::get("CMPIRouting::transferToTarget(messageSize)").print();
372
373        for(int i=0; i<nbTarget; i++)
374        {
375                targetBuffer[i] = new char[targetMessageSize[i]];
376                indexTargetBuffer[i] = 0;
377        }
378
379        char** sourceBuffer = new char*[nbSource];
380        int* indexSourceBuffer = new int[nbSource];
381
382        for(int i=0;i<nbSource; i++)
383        {
384                sourceBuffer[i] = new char[sourceMessageSize[i]];
385                indexSourceBuffer[i]= 0;
386        }
387
388        // pack the data
389        for(int i=0; i<totalTargetElement; i++)
390        {
391                index=targetElementIndex[i];
392                pack(targetElements[i], targetBuffer[index], indexTargetBuffer[index] );
393        }
394
395        indexRequest=0;
396
397        MPI_Barrier(communicator);
398        CTimer::get("CMPIRouting::transferToTarget(data)").reset();
399        CTimer::get("CMPIRouting::transferToTarget(data)").resume();
400        for(int i=0; i<nbSource; i++)
401        {
402                MPI_Irecv(sourceBuffer[i],sourceMessageSize[i],MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
403                indexRequest++;
404        }
405
406        for(int i=0;i<nbTarget; i++)
407        {
408                MPI_Isend(targetBuffer[i],targetMessageSize[i], MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
409                indexRequest++;
410        }
411
412        MPI_Waitall(indexRequest,request,status);
413
414        MPI_Barrier(communicator);
415        CTimer::get("CMPIRouting::transferToTarget(data)").suspend();
416        CTimer::get("CMPIRouting::transferToTarget(data)").print();
417
418        // unpack the data
419        for(int i=0; i<totalSourceElement; i++)
420        {
421                index=sourceElementIndex[i];
422                unpack(sourceElements[i], sourceBuffer[index], indexSourceBuffer[index]);
423        }
424
425        for (int i=0; i<nbTarget; i++) delete[] targetBuffer[i];
426        for (int i=0; i<nbSource; i++) delete[] sourceBuffer[i];
427        delete[] targetBuffer;
428        delete[] indexTargetBuffer;
429        delete[] targetMessageSize;
430        delete[] sourceBuffer;
431        delete[] indexSourceBuffer;
432        delete[] sourceMessageSize;
433        delete[] request;
434        delete[] status;
435}
436
437template<typename T>
438void CMPIRouting::transferFromSource(T* targetElements, T* sourceElements)
439{
440        char** targetBuffer = new char*[nbTarget];
441        int* indexTargetBuffer = new int[nbTarget];
442
443        for (int i = 0; i < nbTarget; i++)
444        {
445                targetBuffer[i] = new char[sizeof(T)*nbTargetElement[i]];
446                indexTargetBuffer[i] = 0;
447        }
448
449        char** sourceBuffer = new char*[nbSource];
450        int* indexSourceBuffer = new int[nbSource];
451
452        for (int i = 0; i < nbSource; i++)
453        {
454                sourceBuffer[i] = new char[sizeof(T)*nbSourceElement[i]];
455                indexSourceBuffer[i] = 0;
456        }
457
458        // pack the data
459        int index;
460        for (int i = 0; i < totalSourceElement; i++)
461        {
462                index = sourceElementIndex[i];
463                *((T*) &(sourceBuffer[index][indexSourceBuffer[index]])) = sourceElements[i];
464                indexSourceBuffer[index] += sizeof(T);
465        }
466
467        MPI_Request* request=new MPI_Request[nbSource+nbTarget];
468        MPI_Status*  status=new MPI_Status[nbSource+nbTarget];
469        int indexRequest=0;
470
471        for(int i=0; i<nbSource; i++)
472        {
473                MPI_Isend(sourceBuffer[i],nbSourceElement[i]*sizeof(T),MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
474                indexRequest++;
475        }
476
477        for(int i=0;i<nbTarget; i++)
478        {
479                MPI_Irecv(targetBuffer[i],nbTargetElement[i]*sizeof(T), MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
480                indexRequest++;
481        }
482
483        MPI_Waitall(indexRequest,request,status);
484
485        // unpack the data
486        for(int i=0;i<totalTargetElement;i++)
487        {
488                index=targetElementIndex[i];
489                targetElements[i] = *((T*) &(targetBuffer[index][indexTargetBuffer[index]]));
490                indexTargetBuffer[index]+=sizeof(T);
491        }
492
493        for(int i=0;i<nbTarget; i++) delete[] targetBuffer[i];
494        for(int i=0;i<nbSource; i++) delete[] sourceBuffer[i];
495        delete[] targetBuffer;
496        delete[] indexTargetBuffer;
497        delete[] sourceBuffer;
498        delete[] indexSourceBuffer;
499        delete[] request;
500        delete[] status;
501}
502
503
504/* number of source and target elements is known from previous call to init() */
505template<typename T, typename t_pack, typename t_unpack>
506void CMPIRouting::transferFromSource(T *targetElements, T *sourceElements, t_pack pack, t_unpack unpack)
507{
508        char **targetBuffer = new char*[nbTarget];
509        int *indexTargetBuffer = new int[nbTarget];
510        int *targetMessageSize = new int[nbTarget];
511        int *sourceMessageSize = new int[nbSource];
512        int index;
513
514        // compute target buffer size
515        for (int i = 0; i < nbSource; i++)  sourceMessageSize[i] = 0;
516
517        for (int i = 0; i < totalSourceElement; i++)
518        {
519                index = sourceElementIndex[i];
520                pack(sourceElements[i], NULL, sourceMessageSize[index]);
521        }
522
523        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
524        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
525        int indexRequest = 0;
526        for (int i = 0; i < nbSource; i++)
527        {
528                MPI_Isend(&sourceMessageSize[i], 1, MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
529                indexRequest++;
530        }
531        for (int i = 0; i < nbTarget; i++)
532        {
533                MPI_Irecv(&targetMessageSize[i], 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
534                indexRequest++;
535        }
536        MPI_Waitall(indexRequest, request, status);
537
538        for (int i = 0; i < nbTarget; i++)
539        {
540                targetBuffer[i] = new char[targetMessageSize[i]];
541                indexTargetBuffer[i] = 0;
542        }
543
544        char **sourceBuffer = new char*[nbSource];
545        int *indexSourceBuffer = new int[nbSource];
546
547        for (int i = 0; i < nbSource; i++)
548        {
549                sourceBuffer[i] = new char[sourceMessageSize[i]];
550                indexSourceBuffer[i] = 0;
551        }
552
553
554        // pack the data
555        for (int i = 0; i < totalSourceElement; i++)
556        {
557                index = sourceElementIndex[i];
558                pack(sourceElements[i], sourceBuffer[index], indexSourceBuffer[index] );
559        }
560
561        indexRequest = 0;
562        for (int i = 0; i < nbSource; i++)
563        {
564                MPI_Isend(sourceBuffer[i], sourceMessageSize[i], MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
565                indexRequest++;
566        }
567        for (int i = 0; i < nbTarget; i++)
568        {
569                MPI_Irecv(targetBuffer[i], targetMessageSize[i], MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
570                indexRequest++;
571        }
572        MPI_Waitall(indexRequest, request, status);
573
574        // unpack the data
575        for (int i = 0; i < totalTargetElement; i++)
576        {
577                index = targetElementIndex[i];
578                unpack(targetElements[i], targetBuffer[index], indexTargetBuffer[index]);
579        }
580
581        for (int i = 0; i < nbTarget; i++) delete[] targetBuffer[i];
582        for (int i = 0; i < nbSource; i++) delete[] sourceBuffer[i];
583        delete[] targetBuffer;
584        delete[] indexTargetBuffer;
585        delete[] targetMessageSize;
586        delete[] sourceBuffer;
587        delete[] indexSourceBuffer;
588        delete[] sourceMessageSize;
589        delete[] request;
590        delete[] status;
591}
592
593CMPIRouting::~CMPIRouting()
594{
595};
596
597template void CMPIRouting::init(const std::vector<int>& route, CMPICascade *cascade);
598template void CMPIRouting::init(const std::vector<vector<int> >& route, CMPICascade *cascade);
599
600template void CMPIRouting::transferToTarget(int *targetElements, int *sourceElements);
601template void CMPIRouting::transferToTarget(int *targetElements, int *sourceElements, void (*pack)(int&, char*, int&), void (* unpack)(int&, char *, int&));
602template void CMPIRouting::transferFromSource(int *targetElements, int *sourceElements);
603template void CMPIRouting::transferFromSource(int *targetElements, int *sourceElements, void (*pack)(int&, char*, int&), void (* unpack)(int&, char *, int&));
604
605template void CMPIRouting::transferToTarget(Node* targetElements, Node* sourceElements, void (*pack)(Node&, char*, int&), void (* unpack)(Node&, char *, int&));
606template void CMPIRouting::transferToTarget(Elt **targetElements, Elt **sourceElements, void (*pack)(Elt *, char*, int&), void (* unpack)(Elt *, char *, int&));
607template void CMPIRouting::transferFromSource(std::vector<int> *targetElements, std::vector<int> *sourceElements, void (*pack)(const std::vector<int>&, char*, int&), void (* unpack)(std::vector<int>&, const char *, int&));
608
609struct NES { int cnt; int risc; int rank; };
610
611template void alltoalls_unknown(const std::vector<std::vector<NES> >& send, std::vector<std::vector<NES> >& recv,
612                                const std::vector<int>& ranks, MPI_Comm communicator);
613
614template void alltoalls_known(const std::vector<std::vector<int> >& send, std::vector<std::vector<int> >& recv,
615                              const std::vector<int>& ranks, MPI_Comm communicator);
616
617}
Note: See TracBrowser for help on using the repository browser.