source: XIOS/dev/dev_trunk_omp/extern/remap/src/mpi_routing.cpp @ 1661

Last change on this file since 1661 was 1661, checked in by yushan, 5 years ago

MARK: branch merged with trunk @1660. Test (test_complete, test_remap) on ADA with IntelMPI and _usingEP/_usingMPI as switch.

File size: 18.6 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
8using namespace ep_lib;
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_Alloc_mem(nbTarget *sizeof(int), MPI_INFO_NULL, &targetRank);
128        MPI_Alloc_mem(nbSource *sizeof(int), MPI_INFO_NULL, &sourceRank);
129
130        targetRankToIndex = new int[mpiSize];
131        int index = 0;
132        for (int i = 0; i < mpiSize; i++)
133        {
134                if (toSend[i] == 1)
135                {
136                        targetRank[index] = i;
137                        targetRankToIndex[i] = index;
138                        index++;
139                }
140        }
141
142        MPI_Barrier(communicator);
143        CTimer::get("CMPIRouting::init(get_source)").reset();
144        CTimer::get("CMPIRouting::init(get_source)").resume();
145
146        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
147        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
148
149        int indexRequest = 0;
150        if (verbose) cout << "CMPIRouting::init     nbSource : " << nbSource << " nbTarget : " << nbTarget << endl;
151//      assert(recvCntDeb == -1 || recvCntDeb == nbSource);
152//cout << mpiRank <<  "DEB : " << recvCntDeb << "    nbSource " << nbSource << " nbTarget : " << nbTarget << endl;
153        for (int i = 0; i < nbSource; i++)
154        {
155                #ifdef _usingEP
156                MPI_Irecv(&sourceRank[i], 1, MPI_INT, -2, 0, communicator, &request[indexRequest++]);
157                #endif
158               
159                #ifdef _usingMPI
160                MPI_Irecv(&sourceRank[i], 1, MPI_INT, MPI_ANY_SOURCE, 0, communicator, &request[indexRequest++]);
161                #endif
162        }
163        MPI_Barrier(communicator);
164        for (int i = 0; i < nbTarget; i++)
165        {
166                MPI_Isend(&mpiRank, 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest++]);
167        }
168        MPI_Waitall(indexRequest, request, status);
169        MPI_Barrier(communicator);  //needed
170        CTimer::get("CMPIRouting::init(get_source)").suspend();
171        CTimer::get("CMPIRouting::init(get_source)").print();
172
173        CTimer::get("CMPIRouting::init(get_source)").reset();
174        CTimer::get("CMPIRouting::init(get_source)").resume();
175
176        indexRequest = 0;
177        for (int i = 0; i < nbSource; i++)
178        {
179#ifdef _usingEP
180                MPI_Irecv(&sourceRank[i], 1, MPI_INT, -2, 0, communicator, &request[indexRequest]);
181#endif
182#ifdef _usingMPI
183                MPI_Irecv(&sourceRank[i], 1, MPI_INT, MPI_ANY_SOURCE, 0, communicator, &request[indexRequest]);
184#endif
185                indexRequest++;
186        }
187
188        for (int i = 0; i < nbTarget; i++)
189        {
190                MPI_Isend(&mpiRank, 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
191                indexRequest++;
192        }
193        MPI_Waitall(indexRequest, request, status);
194        MPI_Barrier(communicator);
195        CTimer::get("CMPIRouting::init(get_source)").suspend();
196        CTimer::get("CMPIRouting::init(get_source)").print();
197
198        CTimer::get("CMPIRouting::init(send_element)").reset();
199        CTimer::get("CMPIRouting::init(send_element)").resume();
200
201        nbTargetElement.resize(nbTarget);
202        nbSourceElement.resize(nbSource);
203
204        for (int i = 0; i < route.size(); i++)
205                setTargetElementIndex(route[i], targetElementIndex, targetRankToIndex);
206
207        for (int i = 0; i < targetElementIndex.size(); i++)
208                nbTargetElement[targetElementIndex[i]]++;
209
210        indexRequest = 0;
211        totalSourceElement = 0;
212        totalTargetElement = 0;
213        for (int i = 0; i < nbSource; i++)
214        {
215                MPI_Irecv(&nbSourceElement[i], 1, MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
216                indexRequest++;
217        }
218
219        for (int i = 0; i < nbTarget; i++)
220        {
221                totalTargetElement += nbTargetElement[i];
222                MPI_Isend(&nbTargetElement[i], 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
223                indexRequest++;
224        }
225
226        MPI_Waitall(indexRequest, request, status);
227
228        CTimer::get("CMPIRouting::init(send_element)").suspend();
229        CTimer::get("CMPIRouting::init(send_element)").print();
230
231        totalSourceElement = 0;
232        for (int i = 0; i < nbSource; i++)
233                totalSourceElement += nbSourceElement[i];
234
235        sourceElementIndex.resize(totalSourceElement);
236
237        totalSourceElement = 0;
238        for (int i = 0; i < nbSource; i++)
239        {
240                for (int j = 0; j < nbSourceElement[i]; j++)
241                {
242                        sourceElementIndex[totalSourceElement] = i;
243                        totalSourceElement++;
244                }
245        }
246
247        delete[]  toSend;
248        delete[]  recvCount;
249        delete[]  request;
250        delete[]  status;
251}
252
253
254int CMPIRouting::getTotalSourceElement(void)
255{
256        return totalSourceElement;
257}
258
259template<typename T>
260void CMPIRouting::transferToTarget(T* targetElements, T* sourceElements)
261{
262        char** targetBuffer = new char*[nbTarget];
263        int* indexTargetBuffer = new int[nbTarget];
264
265        for(int i=0;i<nbTarget; i++)
266        {
267                targetBuffer[i] = new char[sizeof(T)*nbTargetElement[i]];
268                indexTargetBuffer[i]= 0;
269        }
270
271        char** sourceBuffer = new char*[nbSource];
272        int* indexSourceBuffer = new int[nbSource];
273
274        for(int i=0;i<nbSource; i++)
275        {
276                sourceBuffer[i] = new char[sizeof(T)*nbSourceElement[i]];
277                indexSourceBuffer[i]= 0;
278        }
279
280        // pack the data
281        int index;
282        for(int i=0;i<totalTargetElement;i++)
283        {
284                index=targetElementIndex[i];
285                *((T*) &(targetBuffer[index][indexTargetBuffer[index]])) = targetElements[i];
286                indexTargetBuffer[index]+=sizeof(T);
287        }
288
289
290        MPI_Request* request=new MPI_Request[nbSource+nbTarget];
291        MPI_Status*  status=new MPI_Status[nbSource+nbTarget];
292        int indexRequest=0;
293
294        MPI_Barrier(communicator);
295        CTimer::get("CMPIRouting::transferToTarget").reset();
296        CTimer::get("CMPIRouting::transferToTarget").resume();
297
298        for(int i=0; i<nbSource; i++)
299        {
300                MPI_Irecv(sourceBuffer[i],nbSourceElement[i]*sizeof(T),MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
301                indexRequest++;
302        }
303
304        for(int i=0;i<nbTarget; i++)
305        {
306                MPI_Isend(targetBuffer[i],nbTargetElement[i]*sizeof(T), MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
307                indexRequest++;
308        }
309
310        MPI_Waitall(indexRequest,request,status);
311
312        CTimer::get("CMPIRouting::transferToTarget").suspend();
313        CTimer::get("CMPIRouting::transferToTarget").print();
314        MPI_Barrier(communicator);
315
316        // unpack the data
317        for(int i=0;i<totalSourceElement;i++)
318        {
319                index=sourceElementIndex[i];
320                sourceElements[i] = *((T*) &(sourceBuffer[index][indexSourceBuffer[index]]));
321                indexSourceBuffer[index]+=sizeof(T);
322        }
323
324        for(int i=0;i<nbTarget; i++) delete[] targetBuffer[i];
325        for(int i=0;i<nbSource; i++) delete[] sourceBuffer[i];
326        delete[] targetBuffer;
327        delete[] indexTargetBuffer;
328        delete[] sourceBuffer;
329        delete[] indexSourceBuffer;
330        delete[] request;
331        delete[] status;
332}
333
334
335template<typename T, typename t_pack, typename t_unpack>
336void CMPIRouting::transferToTarget(T* targetElements, T* sourceElements, t_pack pack, t_unpack unpack)
337{
338        char** targetBuffer = new char*[nbTarget];
339        int* indexTargetBuffer = new int[nbTarget];
340        int* targetMessageSize = new int[nbTarget];
341        int* sourceMessageSize = new int[nbSource];
342        int index;
343
344        // compute target buffer size
345        for (int i = 0; i < nbTarget; i++)
346                targetMessageSize[i] = 0;
347
348        for (int i = 0; i < totalTargetElement; i++)
349        {
350                index = targetElementIndex[i];
351                pack(targetElements[i], NULL, targetMessageSize[index]);
352        }
353
354        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
355        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
356        int indexRequest = 0;
357
358        MPI_Barrier(communicator);
359        CTimer::get("CMPIRouting::transferToTarget(messageSize)").reset();
360        CTimer::get("CMPIRouting::transferToTarget(messageSize)").resume();
361
362        for(int i=0; i<nbSource; i++)
363        {
364                MPI_Irecv(&sourceMessageSize[i],1,MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
365                indexRequest++;
366        }
367
368        for(int i=0; i<nbTarget; i++)
369        {
370                MPI_Isend(&targetMessageSize[i],1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
371                indexRequest++;
372        }
373
374        MPI_Waitall(indexRequest,request,status);
375
376        MPI_Barrier(communicator);
377        CTimer::get("CMPIRouting::transferToTarget(messageSize)").suspend();
378        CTimer::get("CMPIRouting::transferToTarget(messageSize)").print();
379
380        for(int i=0; i<nbTarget; i++)
381        {
382                targetBuffer[i] = new char[targetMessageSize[i]];
383                indexTargetBuffer[i] = 0;
384        }
385
386        char** sourceBuffer = new char*[nbSource];
387        int* indexSourceBuffer = new int[nbSource];
388
389        for(int i=0;i<nbSource; i++)
390        {
391                sourceBuffer[i] = new char[sourceMessageSize[i]];
392                indexSourceBuffer[i]= 0;
393        }
394
395        // pack the data
396        for(int i=0; i<totalTargetElement; i++)
397        {
398                index=targetElementIndex[i];
399                pack(targetElements[i], targetBuffer[index], indexTargetBuffer[index] );
400        }
401
402        indexRequest=0;
403
404        MPI_Barrier(communicator);
405        CTimer::get("CMPIRouting::transferToTarget(data)").reset();
406        CTimer::get("CMPIRouting::transferToTarget(data)").resume();
407        for(int i=0; i<nbSource; i++)
408        {
409                MPI_Irecv(sourceBuffer[i],sourceMessageSize[i],MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
410                indexRequest++;
411        }
412
413        for(int i=0;i<nbTarget; i++)
414        {
415                MPI_Isend(targetBuffer[i],targetMessageSize[i], MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
416                indexRequest++;
417        }
418
419        MPI_Waitall(indexRequest,request,status);
420
421        MPI_Barrier(communicator);
422        CTimer::get("CMPIRouting::transferToTarget(data)").suspend();
423        CTimer::get("CMPIRouting::transferToTarget(data)").print();
424
425        // unpack the data
426        for(int i=0; i<totalSourceElement; i++)
427        {
428                index=sourceElementIndex[i];
429                unpack(sourceElements[i], sourceBuffer[index], indexSourceBuffer[index]);
430        }
431
432        for (int i=0; i<nbTarget; i++) delete[] targetBuffer[i];
433        for (int i=0; i<nbSource; i++) delete[] sourceBuffer[i];
434        delete[] targetBuffer;
435        delete[] indexTargetBuffer;
436        delete[] targetMessageSize;
437        delete[] sourceBuffer;
438        delete[] indexSourceBuffer;
439        delete[] sourceMessageSize;
440        delete[] request;
441        delete[] status;
442}
443
444template<typename T>
445void CMPIRouting::transferFromSource(T* targetElements, T* sourceElements)
446{
447        char** targetBuffer = new char*[nbTarget];
448        int* indexTargetBuffer = new int[nbTarget];
449
450        for (int i = 0; i < nbTarget; i++)
451        {
452                targetBuffer[i] = new char[sizeof(T)*nbTargetElement[i]];
453                indexTargetBuffer[i] = 0;
454        }
455
456        char** sourceBuffer = new char*[nbSource];
457        int* indexSourceBuffer = new int[nbSource];
458
459        for (int i = 0; i < nbSource; i++)
460        {
461                sourceBuffer[i] = new char[sizeof(T)*nbSourceElement[i]];
462                indexSourceBuffer[i] = 0;
463        }
464
465        // pack the data
466        int index;
467        for (int i = 0; i < totalSourceElement; i++)
468        {
469                index = sourceElementIndex[i];
470                *((T*) &(sourceBuffer[index][indexSourceBuffer[index]])) = sourceElements[i];
471                indexSourceBuffer[index] += sizeof(T);
472        }
473
474        MPI_Request* request=new MPI_Request[nbSource+nbTarget];
475        MPI_Status*  status=new MPI_Status[nbSource+nbTarget];
476        int indexRequest=0;
477
478        for(int i=0; i<nbSource; i++)
479        {
480                MPI_Isend(sourceBuffer[i],nbSourceElement[i]*sizeof(T),MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
481                indexRequest++;
482        }
483
484        for(int i=0;i<nbTarget; i++)
485        {
486                MPI_Irecv(targetBuffer[i],nbTargetElement[i]*sizeof(T), MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
487                indexRequest++;
488        }
489
490        MPI_Waitall(indexRequest,request,status);
491
492        // unpack the data
493        for(int i=0;i<totalTargetElement;i++)
494        {
495                index=targetElementIndex[i];
496                targetElements[i] = *((T*) &(targetBuffer[index][indexTargetBuffer[index]]));
497                indexTargetBuffer[index]+=sizeof(T);
498        }
499
500        for(int i=0;i<nbTarget; i++) delete[] targetBuffer[i];
501        for(int i=0;i<nbSource; i++) delete[] sourceBuffer[i];
502        delete[] targetBuffer;
503        delete[] indexTargetBuffer;
504        delete[] sourceBuffer;
505        delete[] indexSourceBuffer;
506        delete[] request;
507        delete[] status;
508}
509
510
511/* number of source and target elements is known from previous call to init() */
512template<typename T, typename t_pack, typename t_unpack>
513void CMPIRouting::transferFromSource(T *targetElements, T *sourceElements, t_pack pack, t_unpack unpack)
514{
515        char **targetBuffer = new char*[nbTarget];
516        int *indexTargetBuffer = new int[nbTarget];
517        int *targetMessageSize = new int[nbTarget];
518        int *sourceMessageSize = new int[nbSource];
519        int index;
520
521        // compute target buffer size
522        for (int i = 0; i < nbSource; i++)  sourceMessageSize[i] = 0;
523
524        for (int i = 0; i < totalSourceElement; i++)
525        {
526                index = sourceElementIndex[i];
527                pack(sourceElements[i], NULL, sourceMessageSize[index]);
528        }
529
530        MPI_Request *request = new MPI_Request[nbSource + nbTarget];
531        MPI_Status  *status = new MPI_Status[nbSource + nbTarget];
532        int indexRequest = 0;
533        for (int i = 0; i < nbSource; i++)
534        {
535                MPI_Isend(&sourceMessageSize[i], 1, MPI_INT, sourceRank[i], 0, communicator, &request[indexRequest]);
536                indexRequest++;
537        }
538        for (int i = 0; i < nbTarget; i++)
539        {
540                MPI_Irecv(&targetMessageSize[i], 1, MPI_INT, targetRank[i], 0, communicator, &request[indexRequest]);
541                indexRequest++;
542        }
543        MPI_Waitall(indexRequest, request, status);
544
545        for (int i = 0; i < nbTarget; i++)
546        {
547                targetBuffer[i] = new char[targetMessageSize[i]];
548                indexTargetBuffer[i] = 0;
549        }
550
551        char **sourceBuffer = new char*[nbSource];
552        int *indexSourceBuffer = new int[nbSource];
553
554        for (int i = 0; i < nbSource; i++)
555        {
556                sourceBuffer[i] = new char[sourceMessageSize[i]];
557                indexSourceBuffer[i] = 0;
558        }
559
560
561        // pack the data
562        for (int i = 0; i < totalSourceElement; i++)
563        {
564                index = sourceElementIndex[i];
565                pack(sourceElements[i], sourceBuffer[index], indexSourceBuffer[index] );
566        }
567
568        indexRequest = 0;
569        for (int i = 0; i < nbSource; i++)
570        {
571                MPI_Isend(sourceBuffer[i], sourceMessageSize[i], MPI_CHAR, sourceRank[i], 0, communicator, &request[indexRequest]);
572                indexRequest++;
573        }
574        for (int i = 0; i < nbTarget; i++)
575        {
576                MPI_Irecv(targetBuffer[i], targetMessageSize[i], MPI_CHAR, targetRank[i], 0, communicator, &request[indexRequest]);
577                indexRequest++;
578        }
579        MPI_Waitall(indexRequest, request, status);
580
581        // unpack the data
582        for (int i = 0; i < totalTargetElement; i++)
583        {
584                index = targetElementIndex[i];
585                unpack(targetElements[i], targetBuffer[index], indexTargetBuffer[index]);
586        }
587
588        for (int i = 0; i < nbTarget; i++) delete[] targetBuffer[i];
589        for (int i = 0; i < nbSource; i++) delete[] sourceBuffer[i];
590        delete[] targetBuffer;
591        delete[] indexTargetBuffer;
592        delete[] targetMessageSize;
593        delete[] sourceBuffer;
594        delete[] indexSourceBuffer;
595        delete[] sourceMessageSize;
596        delete[] request;
597        delete[] status;
598}
599
600CMPIRouting::~CMPIRouting()
601{
602};
603
604template void CMPIRouting::init(const std::vector<int>& route, CMPICascade *cascade);
605template void CMPIRouting::init(const std::vector<vector<int> >& route, CMPICascade *cascade);
606
607template void CMPIRouting::transferToTarget(int *targetElements, int *sourceElements);
608template void CMPIRouting::transferToTarget(int *targetElements, int *sourceElements, void (*pack)(int&, char*, int&), void (* unpack)(int&, char *, int&));
609template void CMPIRouting::transferFromSource(int *targetElements, int *sourceElements);
610template void CMPIRouting::transferFromSource(int *targetElements, int *sourceElements, void (*pack)(int&, char*, int&), void (* unpack)(int&, char *, int&));
611
612template void CMPIRouting::transferToTarget(Node* targetElements, Node* sourceElements, void (*pack)(Node&, char*, int&), void (* unpack)(Node&, char *, int&));
613template void CMPIRouting::transferToTarget(Elt **targetElements, Elt **sourceElements, void (*pack)(Elt *, char*, int&), void (* unpack)(Elt *, char *, int&));
614template 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&));
615
616struct NES { int cnt; int risc; int rank; };
617
618template void alltoalls_unknown(const std::vector<std::vector<NES> >& send, std::vector<std::vector<NES> >& recv,
619                                const std::vector<int>& ranks, MPI_Comm communicator);
620
621template void alltoalls_known(const std::vector<std::vector<int> >& send, std::vector<std::vector<int> >& recv,
622                              const std::vector<int>& ranks, MPI_Comm communicator);
623
624}
Note: See TracBrowser for help on using the repository browser.