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