source: XIOS/dev/branch_openmp/extern/remap/src/mpi_routing.cpp @ 1362

Last change on this file since 1362 was 1362, checked in by yushan, 6 years ago

unify type : MPI_Message MPI_Info

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