source: XIOS/trunk/extern/remap/src/mapper.cpp @ 815

Last change on this file since 815 was 694, checked in by rlacroix, 9 years ago

Fix compilation issues caused by the new "remap" library.

Use our MPI header so that C++ bindings are always disabled.

File size: 27.7 KB
Line 
1#include "mpi.hpp"
2#include <map>
3#include "cputime.hpp"       /* time */
4#include "meshutil.hpp"
5#include "polyg.hpp"
6#include "circle.hpp"
7#include "intersect.hpp"
8#include "intersection_ym.hpp"
9#include "errhandle.hpp"
10#include "mpi_routing.hpp"
11#include "gridRemap.hpp"
12
13#include "mapper.hpp"
14
15namespace sphereRemap {
16
17/* A subdivition of an array into N sub-arays
18   can be represented by the length of the N arrays
19   or by the offsets when each of the N arrays starts.
20   This function convertes from the former to the latter. */
21void cptOffsetsFromLengths(const int *lengths, int *offsets, int sz)
22{
23        offsets[0] = 0;
24        for (int i = 1; i < sz; i++)
25                offsets[i] = offsets[i-1] + lengths[i-1];
26}
27
28
29void Mapper::setSourceMesh(const double* boundsLon, const double* boundsLat, int nVertex, int nbCells, const double* pole, const long int* globalId)
30{
31  srcGrid.pole = Coord(pole[0], pole[1], pole[2]);
32
33        int mpiRank, mpiSize;
34        MPI_Comm_rank(communicator, &mpiRank);
35        MPI_Comm_size(communicator, &mpiSize);
36
37        sourceElements.reserve(nbCells);
38        sourceMesh.reserve(nbCells);
39  sourceGlobalId.resize(nbCells) ;
40
41  if (globalId==NULL)
42  {
43    long int offset ;
44    long int nb=nbCells ;
45    MPI_Scan(&nb,&offset,1,MPI_LONG,MPI_SUM,communicator) ;
46    offset=offset-nb ;
47    for(int i=0;i<nbCells;i++) sourceGlobalId[i]=offset+i ;
48  }
49  else sourceGlobalId.assign(globalId,globalId+nbCells);
50
51        for (int i = 0; i < nbCells; i++)
52        {
53                int offs = i*nVertex;
54                Elt elt(boundsLon + offs, boundsLat + offs, nVertex);
55                elt.src_id.rank = mpiRank;
56                elt.src_id.ind = i;
57    elt.src_id.globalId = sourceGlobalId[i];
58                sourceElements.push_back(elt);
59                sourceMesh.push_back(Node(elt.x, cptRadius(elt), &sourceElements.back()));
60                cptEltGeom(sourceElements[i], Coord(pole[0], pole[1], pole[2]));
61        }
62
63}
64
65void Mapper::setTargetMesh(const double* boundsLon, const double* boundsLat, int nVertex, int nbCells, const double* pole, const long int* globalId)
66{
67  tgtGrid.pole = Coord(pole[0], pole[1], pole[2]);
68
69        int mpiRank, mpiSize;
70        MPI_Comm_rank(communicator, &mpiRank);
71        MPI_Comm_size(communicator, &mpiSize);
72
73        targetElements.reserve(nbCells);
74        targetMesh.reserve(nbCells);
75
76  targetGlobalId.resize(nbCells) ;
77  if (globalId==NULL)
78  {
79    long int offset ;
80    long int nb=nbCells ;
81    MPI_Scan(&nb,&offset,1,MPI_LONG,MPI_SUM,communicator) ;
82    offset=offset-nb ;
83    for(int i=0;i<nbCells;i++) targetGlobalId[i]=offset+i ;
84  }
85  else targetGlobalId.assign(globalId,globalId+nbCells);
86
87        for (int i = 0; i < nbCells; i++)
88        {
89                int offs = i*nVertex;
90                Elt elt(boundsLon + offs, boundsLat + offs, nVertex);
91                targetElements.push_back(elt);
92                targetMesh.push_back(Node(elt.x, cptRadius(elt), &sourceElements.back()));
93                cptEltGeom(targetElements[i], Coord(pole[0], pole[1], pole[2]));
94        }
95
96
97}
98
99void Mapper::setSourceValue(const double* val)
100{
101  int size=sourceElements.size() ;
102  for(int i=0;i<size;++i) sourceElements[i].val=val[i] ;
103}
104
105void Mapper::getTargetValue(double* val)
106{
107  int size=targetElements.size() ;
108  for(int i=0;i<size;++i) val[i]=targetElements[i].val ;
109}
110
111vector<double> Mapper::computeWeights(int interpOrder)
112{
113        vector<double> timings;
114        int mpiSize, mpiRank;
115        MPI_Comm_size(communicator, &mpiSize);
116        MPI_Comm_rank(communicator, &mpiRank);
117
118  this->buildSSTree(sourceMesh, targetMesh);
119
120        if (mpiRank == 0 && verbose) cout << "Computing intersections ..." << endl;
121        double tic = cputime();
122        computeIntersection(&targetElements[0], targetElements.size());
123        timings.push_back(cputime() - tic);
124
125        tic = cputime();
126        if (interpOrder == 2) {
127                if (mpiRank == 0 && verbose) cout << "Computing grads ..." << endl;
128                buildMeshTopology();
129                computeGrads();
130        }
131        timings.push_back(cputime() - tic);
132
133        /* Prepare computation of weights */
134        /* compute number of intersections which for the first order case
135           corresponds to the number of edges in the remap matrix */
136        int nIntersections = 0;
137        for (int j = 0; j < targetElements.size(); j++)
138        {
139                Elt &elt = targetElements[j];
140                for (list<Polyg*>::iterator it = elt.is.begin(); it != elt.is.end(); it++)
141                        nIntersections++;
142        }
143        /* overallocate for NMAX neighbours for each elements */
144        remapMatrix = new double[nIntersections*NMAX];
145        srcAddress = new int[nIntersections*NMAX];
146        srcRank = new int[nIntersections*NMAX];
147        dstAddress = new int[nIntersections*NMAX];
148  sourceWeightId =new long[nIntersections*NMAX];
149  targetWeightId =new long[nIntersections*NMAX];
150
151
152        if (mpiRank == 0 && verbose) cout << "Remapping..." << endl;
153        tic = cputime();
154        nWeights = remap(&targetElements[0], targetElements.size(), interpOrder);
155        timings.push_back(cputime() - tic);
156
157  for (int i = 0; i < targetElements.size(); i++) targetElements[i].delete_intersections();
158
159        return timings;
160}
161
162/**
163   @param elements are cells of the target grid that are distributed over CPUs
164          indepentently of the distribution of the SS-tree.
165   @param nbElements is the size of the elements array.
166   @param order is the order of interpolaton (must be 1 or 2).
167*/
168int Mapper::remap(Elt *elements, int nbElements, int order)
169{
170        int mpiSize, mpiRank;
171        MPI_Comm_size(communicator, &mpiSize);
172        MPI_Comm_rank(communicator, &mpiRank);
173
174        /* create list of intersections (super mesh elements) for each rank */
175        multimap<int, Polyg *> *elementList = new multimap<int, Polyg *>[mpiSize];
176        for (int j = 0; j < nbElements; j++)
177        {
178                Elt& e = elements[j];
179                for (list<Polyg *>::iterator it = e.is.begin(); it != e.is.end(); it++)
180                        elementList[(*it)->id.rank].insert(pair<int, Polyg *>((*it)->id.ind, *it));
181        }
182
183        int *nbSendElement = new int[mpiSize];
184        int **sendElement = new int*[mpiSize]; /* indices of elements required from other rank */
185        double **recvValue = new double*[mpiSize];
186        Coord **recvGrad = new Coord*[mpiSize];
187        GloId **recvNeighIds = new GloId*[mpiSize]; /* ids of the of the source neighbours which also contribute through gradient */
188        for (int rank = 0; rank < mpiSize; rank++)
189        {
190                /* get size for allocation */
191                int last = -1; /* compares unequal to any index */
192                int index = -1; /* increased to starting index 0 in first iteration */
193                for (multimap<int, Polyg *>::iterator it = elementList[rank].begin(); it != elementList[rank].end(); ++it)
194                {
195                        if (last != it->first)
196                                index++;
197                        (it->second)->id.ind = index;
198                        last = it->first;
199                }
200                nbSendElement[rank] = index + 1;
201
202                /* if size is non-zero allocate and collect indices of elements on other ranks that we intersect */
203                if (nbSendElement[rank] > 0)
204                {
205                        sendElement[rank] = new int[nbSendElement[rank]];
206                        recvValue[rank]   = new double[nbSendElement[rank]];
207                        if (order == 2)
208                        {
209                                recvNeighIds[rank] = new GloId[nbSendElement[rank]*(NMAX+1)];
210                                recvGrad[rank]    = new Coord[nbSendElement[rank]*(NMAX+1)];
211                        }
212                        else
213                                recvNeighIds[rank] = new GloId[nbSendElement[rank]];
214
215                        last = -1;
216                        index = -1;
217                        for (multimap<int, Polyg *>::iterator it = elementList[rank].begin(); it != elementList[rank].end(); ++it)
218                        {
219                                if (last != it->first)
220                                        index++;
221                                sendElement[rank][index] = it->first;
222                                last = it->first;
223                        }
224                }
225        }
226
227        /* communicate sizes of source elements to be sent (index lists and later values and gradients) */
228        int *nbRecvElement = new int[mpiSize];
229        MPI_Alltoall(nbSendElement, 1, MPI_INT, nbRecvElement, 1, MPI_INT, communicator);
230
231        /* communicate indices of source elements on other ranks whoes value and gradient we need (since intersection) */
232        int nbSendRequest = 0;
233        int nbRecvRequest = 0;
234        int **recvElement = new int*[mpiSize];
235        double **sendValue = new double*[mpiSize];
236        Coord **sendGrad = new Coord*[mpiSize];
237        GloId **sendNeighIds = new GloId*[mpiSize];
238        MPI_Request *sendRequest = new MPI_Request[3*mpiSize];
239        MPI_Request *recvRequest = new MPI_Request[3*mpiSize];
240        for (int rank = 0; rank < mpiSize; rank++)
241        {
242                if (nbSendElement[rank] > 0)
243                {
244                        MPI_Issend(sendElement[rank], nbSendElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
245                        nbSendRequest++;
246                }
247
248                if (nbRecvElement[rank] > 0)
249                {
250                        recvElement[rank] = new int[nbRecvElement[rank]];
251                        sendValue[rank]   = new double[nbRecvElement[rank]];
252                        if (order == 2)
253                        {
254                                sendNeighIds[rank] = new GloId[nbRecvElement[rank]*(NMAX+1)];
255                                sendGrad[rank]    = new Coord[nbRecvElement[rank]*(NMAX+1)];
256                        }
257                        else
258                        {
259                                sendNeighIds[rank] = new GloId[nbRecvElement[rank]];
260                        }
261                        MPI_Irecv(recvElement[rank], nbRecvElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
262                        nbRecvRequest++;
263                }
264        }
265        MPI_Status *status = new MPI_Status[3*mpiSize];
266        MPI_Waitall(nbRecvRequest, recvRequest, status);
267        MPI_Waitall(nbSendRequest, sendRequest, status);
268
269        /* for all indices that have been received from requesting ranks: pack values and gradients, then send */
270        nbSendRequest = 0;
271        nbRecvRequest = 0;
272        for (int rank = 0; rank < mpiSize; rank++)
273        {
274                if (nbRecvElement[rank] > 0)
275                {
276                        int jj = 0; // jj == j if no weight writing
277                        for (int j = 0; j < nbRecvElement[rank]; j++)
278                        {
279                                sendValue[rank][j] = sstree.localElements[recvElement[rank][j]].val;
280                                if (order == 2)
281                                {
282                                        sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].grad;
283                                        sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].src_id;
284                                        jj++;
285                                        for (int i = 0; i < NMAX; i++)
286                                        {
287                                                sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].gradNeigh[i];
288                                                sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].neighId[i];
289                                                jj++;
290                                        }
291                                }
292                                else
293                                        sendNeighIds[rank][j] = sstree.localElements[recvElement[rank][j]].src_id;
294                        }
295                        MPI_Issend(sendValue[rank],  nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]);
296                        nbSendRequest++;
297                        if (order == 2)
298                        {
299                                MPI_Issend(sendGrad[rank], 3*nbRecvElement[rank]*(NMAX+1),
300                                                                MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]);
301                                nbSendRequest++;
302                                MPI_Issend(sendNeighIds[rank], 4*nbRecvElement[rank]*(NMAX+1), MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
303//ym  --> attention taille GloId
304                                nbSendRequest++;
305                        }
306                        else
307                        {
308                                MPI_Issend(sendNeighIds[rank], 4*nbRecvElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
309//ym  --> attention taille GloId
310                                nbSendRequest++;
311                        }
312                }
313                if (nbSendElement[rank] > 0)
314                {
315                        MPI_Irecv(recvValue[rank],  nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]);
316                        nbRecvRequest++;
317                        if (order == 2)
318                        {
319                                MPI_Irecv(recvGrad[rank], 3*nbSendElement[rank]*(NMAX+1),
320                                                MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]);
321                                nbRecvRequest++;
322                                MPI_Irecv(recvNeighIds[rank], 4*nbSendElement[rank]*(NMAX+1), MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
323//ym  --> attention taille GloId
324                                nbRecvRequest++;
325                        }
326                        else
327                        {
328                                MPI_Irecv(recvNeighIds[rank], 4*nbSendElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
329//ym  --> attention taille GloId
330                                nbRecvRequest++;
331                        }
332                }
333        }
334        MPI_Waitall(nbRecvRequest, recvRequest, status);
335        MPI_Waitall(nbSendRequest, sendRequest, status);
336
337        /* now that all values and gradients are available use them to computed interpolated values on target
338           and also to compute weights */
339        int i = 0;
340        for (int j = 0; j < nbElements; j++)
341        {
342                Elt& e = elements[j];
343
344                /* since for the 2nd order case source grid elements can contribute to a destination grid element over several "paths"
345                   (step1: gradient is computed using neighbours on same grid, step2: intersection uses several elements on other grid)
346                   accumulate them so that there is only one final weight between two elements */
347                map<GloId,double> wgt_map;
348
349                /* for destination element `e` loop over all intersetions/the corresponding source elements */
350                for (list<Polyg *>::iterator it = e.is.begin(); it != e.is.end(); it++)
351                {
352                        /* it is the intersection element, so it->x and it->area are barycentre and area of intersection element (super mesh)
353                        but it->id is id of the source element that it intersects */
354                        int n1 = (*it)->id.ind;
355                        int rank = (*it)->id.rank;
356                        double fk = recvValue[rank][n1];
357                        double w = (*it)->area;
358
359                        /* first order: src value times weight (weight = supermesh area), later divide by target area */
360                        int kk = (order == 2) ? n1 * (NMAX + 1) : n1;
361                        GloId neighID = recvNeighIds[rank][kk];
362                        wgt_map[neighID] += (*it)->area;
363
364                        if (order == 2)
365                        {
366                                for (int k = 0; k < NMAX+1; k++)
367                                {
368                                        int kk = n1 * (NMAX + 1) + k;
369                                        GloId neighID = recvNeighIds[rank][kk];
370                                        if (neighID.ind == -1) break;
371                                        wgt_map[neighID] +=
372                                                w * scalarprod(recvGrad[rank][kk], (*it)->x);
373                                }
374
375                        }
376                }
377                for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++)
378                {
379                        this->remapMatrix[i] = it->second / e.area;
380                        this->srcAddress[i] = it->first.ind;
381                        this->srcRank[i] = it->first.rank;
382                        this->dstAddress[i] = j;
383      this->sourceWeightId[i]= it->first.globalId ;
384      this->targetWeightId[i]= targetGlobalId[j] ;
385                        i++;
386                }
387        }
388
389        /* free all memory allocated in this function */
390        for (int rank = 0; rank < mpiSize; rank++)
391        {
392                if (nbSendElement[rank] > 0)
393                {
394                        delete[] sendElement[rank];
395                        delete[] recvValue[rank];
396                        if (order == 2)
397                        {
398                                delete[] recvGrad[rank];
399                        }
400                        delete[] recvNeighIds[rank];
401                }
402                if (nbRecvElement[rank] > 0)
403                {
404                        delete[] recvElement[rank];
405                        delete[] sendValue[rank];
406                        if (order == 2)
407                                delete[] sendGrad[rank];
408                        delete[] sendNeighIds[rank];
409                }
410        }
411        delete[] status;
412        delete[] sendRequest;
413        delete[] recvRequest;
414        delete[] elementList;
415        delete[] nbSendElement;
416        delete[] nbRecvElement;
417        delete[] sendElement;
418        delete[] recvElement;
419        delete[] sendValue;
420        delete[] recvValue;
421        delete[] sendGrad;
422        delete[] recvGrad;
423        delete[] sendNeighIds;
424        delete[] recvNeighIds;
425        return i;
426}
427
428void Mapper::computeGrads()
429{
430        /* array of pointers to collect local elements and elements received from other cpu */
431        vector<Elt*> globalElements(sstree.nbLocalElements + nbNeighbourElements);
432        int index = 0;
433        for (int i = 0; i < sstree.nbLocalElements; i++, index++)
434                globalElements[index] = &(sstree.localElements[i]);
435        for (int i = 0; i < nbNeighbourElements; i++, index++)
436                globalElements[index] = &neighbourElements[i];
437
438        update_baryc(sstree.localElements, sstree.nbLocalElements);
439        computeGradients(&globalElements[0], sstree.nbLocalElements);
440}
441
442/** for each element of the source grid, finds all the neighbouring elements that share an edge
443    (filling array neighbourElements). This is used later to compute gradients */
444void Mapper::buildMeshTopology()
445{
446        int mpiSize, mpiRank;
447        MPI_Comm_size(communicator, &mpiSize);
448        MPI_Comm_rank(communicator, &mpiRank);
449
450        vector<Node> *routingList = new vector<Node>[mpiSize];
451        vector<vector<int> > routes(sstree.localTree.leafs.size());
452
453        sstree.routeIntersections(routes, sstree.localTree.leafs);
454
455        for (int i = 0; i < routes.size(); ++i)
456                for (int k = 0; k < routes[i].size(); ++k)
457                        routingList[routes[i][k]].push_back(sstree.localTree.leafs[i]);
458        routingList[mpiRank].clear();
459
460
461        CMPIRouting mpiRoute(communicator);
462        mpiRoute.init(routes);
463        int nRecv = mpiRoute.getTotalSourceElement();
464        cout << mpiRank << " NRECV " << nRecv << "(" << routes.size() << ")"<< endl;
465
466        int *nbSendNode = new int[mpiSize];
467        int *nbRecvNode = new int[mpiSize];
468        int *sendMessageSize = new int[mpiSize];
469        int *recvMessageSize = new int[mpiSize];
470
471        for (int rank = 0; rank < mpiSize; rank++)
472        {
473                nbSendNode[rank] = routingList[rank].size();
474                sendMessageSize[rank] = 0;
475                for (size_t j = 0; j < routingList[rank].size(); j++)
476                {
477                        Elt *elt = (Elt *) (routingList[rank][j].data);
478                        sendMessageSize[rank] += packedPolygonSize(*elt);
479                }
480        }
481
482        MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
483        MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
484
485        char **sendBuffer = new char*[mpiSize];
486        char **recvBuffer = new char*[mpiSize];
487        int *pos = new int[mpiSize];
488
489        for (int rank = 0; rank < mpiSize; rank++)
490        {
491                if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sendMessageSize[rank]];
492                if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]];
493        }
494
495        for (int rank = 0; rank < mpiSize; rank++)
496        {
497                pos[rank] = 0;
498                for (size_t j = 0; j < routingList[rank].size(); j++)
499                {
500                        Elt *elt = (Elt *) (routingList[rank][j].data);
501                        packPolygon(*elt, sendBuffer[rank], pos[rank]);
502                }
503        }
504        delete [] routingList;
505
506
507        int nbSendRequest = 0;
508        int nbRecvRequest = 0;
509        MPI_Request *sendRequest = new MPI_Request[mpiSize];
510        MPI_Request *recvRequest = new MPI_Request[mpiSize];
511        MPI_Status  *status      = new MPI_Status[mpiSize];
512
513        for (int rank = 0; rank < mpiSize; rank++)
514        {
515                if (nbSendNode[rank] > 0)
516                {
517                        MPI_Issend(sendBuffer[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
518                        nbSendRequest++;
519                }
520                if (nbRecvNode[rank] > 0)
521                {
522                        MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
523                        nbRecvRequest++;
524                }
525        }
526
527        MPI_Waitall(nbRecvRequest, recvRequest, status);
528        MPI_Waitall(nbSendRequest, sendRequest, status);
529
530        for (int rank = 0; rank < mpiSize; rank++)
531                if (nbSendNode[rank] > 0) delete [] sendBuffer[rank];
532        delete [] sendBuffer;
533
534        char **sendBuffer2 = new char*[mpiSize];
535        char **recvBuffer2 = new char*[mpiSize];
536
537        for (int rank = 0; rank < mpiSize; rank++)
538        {
539                nbSendNode[rank] = 0;
540                sendMessageSize[rank] = 0;
541
542                if (nbRecvNode[rank] > 0)
543                {
544                        set<NodePtr> neighbourList;
545                        pos[rank] = 0;
546                        for (int j = 0; j < nbRecvNode[rank]; j++)
547                        {
548                                Elt elt;
549                                unpackPolygon(elt, recvBuffer[rank], pos[rank]);
550                                Node node(elt.x, cptRadius(elt), &elt);
551                                findNeighbour(sstree.localTree.root, &node, neighbourList);
552                        }
553                        nbSendNode[rank] = neighbourList.size();
554                        for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++)
555                        {
556                                Elt *elt = (Elt *) ((*it)->data);
557                                sendMessageSize[rank] += packedPolygonSize(*elt);
558                        }
559
560                        sendBuffer2[rank] = new char[sendMessageSize[rank]];
561                        pos[rank] = 0;
562
563                        for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++)
564                        {
565                                Elt *elt = (Elt *) ((*it)->data);
566                                packPolygon(*elt, sendBuffer2[rank], pos[rank]);
567                        }
568                }
569        }
570        for (int rank = 0; rank < mpiSize; rank++)
571                if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank];
572        delete [] recvBuffer;
573
574
575        MPI_Barrier(communicator);
576        MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
577        MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
578
579        for (int rank = 0; rank < mpiSize; rank++)
580                if (nbRecvNode[rank] > 0) recvBuffer2[rank] = new char[recvMessageSize[rank]];
581
582        nbSendRequest = 0;
583        nbRecvRequest = 0;
584
585        for (int rank = 0; rank < mpiSize; rank++)
586        {
587                if (nbSendNode[rank] > 0)
588                {
589                        MPI_Issend(sendBuffer2[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
590                        nbSendRequest++;
591                }
592                if (nbRecvNode[rank] > 0)
593                {
594                        MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
595                        nbRecvRequest++;
596                }
597        }
598
599        MPI_Waitall(nbRecvRequest, recvRequest, status);
600        MPI_Waitall(nbSendRequest, sendRequest, status);
601
602        int nbNeighbourNodes = 0;
603        for (int rank = 0; rank < mpiSize; rank++)
604                nbNeighbourNodes += nbRecvNode[rank];
605
606        neighbourElements = new Elt[nbNeighbourNodes];
607        nbNeighbourElements = nbNeighbourNodes;
608
609        int index = 0;
610        for (int rank = 0; rank < mpiSize; rank++)
611        {
612                pos[rank] = 0;
613                for (int j = 0; j < nbRecvNode[rank]; j++)
614                {
615                        unpackPolygon(neighbourElements[index], recvBuffer2[rank], pos[rank]);
616                        neighbourElements[index].id.ind = sstree.localTree.leafs.size() + index;
617                        index++;
618                }
619        }
620        for (int rank = 0; rank < mpiSize; rank++)
621        {
622                if (nbRecvNode[rank] > 0) delete [] recvBuffer2[rank];
623                if (nbSendNode[rank] > 0) delete [] sendBuffer2[rank];
624        }
625        delete [] recvBuffer2;
626        delete [] sendBuffer2;
627        delete [] sendMessageSize;
628        delete [] recvMessageSize;
629        delete [] nbSendNode;
630        delete [] nbRecvNode;
631        delete [] sendRequest;
632        delete [] recvRequest;
633        delete [] status;
634        delete [] pos;
635
636        /* re-compute on received elements to avoid having to send this information */
637        neighbourNodes.resize(nbNeighbourNodes);
638        setCirclesAndLinks(neighbourElements, neighbourNodes);
639        cptAllEltsGeom(neighbourElements, nbNeighbourNodes, srcGrid.pole);
640
641        /* the local SS tree must include nodes from other cpus if they are potential
642           intersector of a local node */
643        sstree.localTree.insertNodes(neighbourNodes);
644
645        /* for every local element,
646           use the SS-tree to find all elements (including neighbourElements)
647           who are potential neighbours because their circles intersect,
648           then check all canditates for common edges to build up connectivity information
649        */
650        for (int j = 0; j < sstree.localTree.leafs.size(); j++)
651        {
652                Node& node = sstree.localTree.leafs[j];
653
654                /* find all leafs whoes circles that intersect node's circle and save into node->intersectors */
655                node.search(sstree.localTree.root);
656
657                Elt *elt = (Elt *)(node.data);
658
659                for (int i = 0; i < elt->n; i++) elt->neighbour[i] = NOT_FOUND;
660
661                /* for element `elt` loop through all nodes in the SS-tree
662                   whoes circles intersect with the circle around `elt` (the SS intersectors)
663                   and check if they are neighbours in the sense that the two elements share an edge.
664                   If they do, save this information for elt */
665                for (list<NodePtr>::iterator it = (node.intersectors).begin(); it != (node.intersectors).end(); ++it)
666                {
667                        Elt *elt2 = (Elt *)((*it)->data);
668                        set_neighbour(*elt, *elt2);
669                }
670
671                for (int i = 0; i < elt->n; i++)
672                {
673                        if (elt->neighbour[i] == NOT_FOUND)
674                                error_exit("neighbour not found");
675                }
676        }
677}
678
679/** @param elements are the target grid elements */
680void Mapper::computeIntersection(Elt *elements, int nbElements)
681{
682        int mpiSize, mpiRank;
683        MPI_Comm_size(communicator, &mpiSize);
684        MPI_Comm_rank(communicator, &mpiRank);
685
686        MPI_Barrier(communicator);
687
688        vector<Node> *routingList = new vector<Node>[mpiSize];
689
690        vector<Node> routeNodes;  routeNodes.reserve(nbElements);
691        for (int j = 0; j < nbElements; j++)
692        {
693                elements[j].id.ind = j;
694                elements[j].id.rank = mpiRank;
695                routeNodes.push_back(Node(elements[j].x, cptRadius(elements[j]), &elements[j]));
696        }
697
698        vector<vector<int> > routes(routeNodes.size());
699        sstree.routeIntersections(routes, routeNodes);
700        for (int i = 0; i < routes.size(); ++i)
701                for (int k = 0; k < routes[i].size(); ++k)
702                        routingList[routes[i][k]].push_back(routeNodes[i]);
703
704        if (verbose >= 2)
705        {
706                cout << " --> rank  " << mpiRank << " nbElements " << nbElements << " : ";
707                for (int rank = 0; rank < mpiSize; rank++)
708                        cout << routingList[rank].size() << "   ";
709                cout << endl;
710        }
711        MPI_Barrier(communicator);
712
713        int *nbSendNode = new int[mpiSize];
714        int *nbRecvNode = new int[mpiSize];
715        int *sentMessageSize = new int[mpiSize];
716        int *recvMessageSize = new int[mpiSize];
717
718        for (int rank = 0; rank < mpiSize; rank++)
719        {
720                nbSendNode[rank] = routingList[rank].size();
721                sentMessageSize[rank] = 0;
722                for (size_t j = 0; j < routingList[rank].size(); j++)
723                {
724                        Elt *elt = (Elt *) (routingList[rank][j].data);
725                        sentMessageSize[rank] += packedPolygonSize(*elt);
726                }
727        }
728
729        MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
730        MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
731
732        int total = 0;
733
734        for (int rank = 0; rank < mpiSize; rank++)
735        {
736                total = total + nbRecvNode[rank];
737        }
738
739        if (verbose >= 2) cout << "---> rank " << mpiRank << " : compute intersection : total received nodes  " << total << endl;
740        char **sendBuffer = new char*[mpiSize];
741        char **recvBuffer = new char*[mpiSize];
742        int *pos = new int[mpiSize];
743
744        for (int rank = 0; rank < mpiSize; rank++)
745        {
746                if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sentMessageSize[rank]];
747                if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]];
748        }
749
750        for (int rank = 0; rank < mpiSize; rank++)
751        {
752                pos[rank] = 0;
753                for (size_t j = 0; j < routingList[rank].size(); j++)
754                {
755                        Elt* elt = (Elt *) (routingList[rank][j].data);
756                        packPolygon(*elt, sendBuffer[rank], pos[rank]);
757                }
758        }
759        delete [] routingList;
760
761        int nbSendRequest = 0;
762        int nbRecvRequest = 0;
763        MPI_Request *sendRequest = new MPI_Request[mpiSize];
764        MPI_Request *recvRequest = new MPI_Request[mpiSize];
765        MPI_Status   *status = new MPI_Status[mpiSize];
766
767        for (int rank = 0; rank < mpiSize; rank++)
768        {
769                if (nbSendNode[rank] > 0)
770                {
771                        MPI_Issend(sendBuffer[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
772                        nbSendRequest++;
773                }
774                if (nbRecvNode[rank] > 0)
775                {
776                        MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
777                        nbRecvRequest++;
778                }
779        }
780
781        MPI_Waitall(nbRecvRequest, recvRequest, status);
782        MPI_Waitall(nbSendRequest, sendRequest, status);
783        char **sendBuffer2 = new char*[mpiSize];
784        char **recvBuffer2 = new char*[mpiSize];
785
786        double tic = cputime();
787        for (int rank = 0; rank < mpiSize; rank++)
788        {
789                sentMessageSize[rank] = 0;
790
791                if (nbRecvNode[rank] > 0)
792                {
793                        Elt *recvElt = new Elt[nbRecvNode[rank]];
794                        pos[rank] = 0;
795                        for (int j = 0; j < nbRecvNode[rank]; j++)
796                        {
797                                unpackPolygon(recvElt[j], recvBuffer[rank], pos[rank]);
798                                cptEltGeom(recvElt[j], tgtGrid.pole);
799                                Node recvNode(recvElt[j].x, cptRadius(recvElt[j]), &recvElt[j]);
800                                recvNode.search(sstree.localTree.root);
801                                /* for a node holding an element of the target, loop throught candidates for intersecting source */
802                                for (list<NodePtr>::iterator it = (recvNode.intersectors).begin(); it != (recvNode.intersectors).end(); ++it)
803                                {
804                                        Elt *elt2 = (Elt *) ((*it)->data);
805                                        /* recvElt is target, elt2 is source */
806//                                      intersect(&recvElt[j], elt2);
807                                        intersect_ym(&recvElt[j], elt2);
808                                }
809
810                                if (recvElt[j].is.size() > 0) sentMessageSize[rank] += packIntersectionSize(recvElt[j]);
811
812                                // here recvNode goes out of scope
813                        }
814
815                        if (sentMessageSize[rank] > 0)
816                        {
817                                sentMessageSize[rank] += sizeof(int);
818                                sendBuffer2[rank] = new char[sentMessageSize[rank]];
819                                *((int *) sendBuffer2[rank]) = 0;
820                                pos[rank] = sizeof(int);
821                                for (int j = 0; j < nbRecvNode[rank]; j++)
822                                {
823                                        packIntersection(recvElt[j], sendBuffer2[rank], pos[rank]);
824                                        //FIXME should be deleted: recvElt[j].delete_intersections(); // intersection areas have been packed to buffer and won't be used any more
825                                }
826                        }
827                        delete [] recvElt;
828
829                }
830        }
831        delete [] pos;
832
833        for (int rank = 0; rank < mpiSize; rank++)
834        {
835                if (nbSendNode[rank] > 0) delete [] sendBuffer[rank];
836                if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank];
837                nbSendNode[rank] = 0;
838        }
839
840        if (verbose >= 2) cout << "Rank " << mpiRank << "  Compute (internal) intersection " << cputime() - tic << " s" << endl;
841        MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
842
843        for (int rank = 0; rank < mpiSize; rank++)
844                if (recvMessageSize[rank] > 0)
845                        recvBuffer2[rank] = new char[recvMessageSize[rank]];
846
847        nbSendRequest = 0;
848        nbRecvRequest = 0;
849
850        for (int rank = 0; rank < mpiSize; rank++)
851        {
852                if (sentMessageSize[rank] > 0)
853                {
854                        MPI_Issend(sendBuffer2[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
855                        nbSendRequest++;
856                }
857                if (recvMessageSize[rank] > 0)
858                {
859                        MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
860                        nbRecvRequest++;
861                }
862        }
863
864        MPI_Waitall(nbRecvRequest, recvRequest, status);
865        MPI_Waitall(nbSendRequest, sendRequest, status);
866
867        delete [] sendRequest;
868        delete [] recvRequest;
869        delete [] status;
870        for (int rank = 0; rank < mpiSize; rank++)
871        {
872                if (nbRecvNode[rank] > 0)
873                {
874                        if (sentMessageSize[rank] > 0)
875                                delete [] sendBuffer2[rank];
876                }
877
878                if (recvMessageSize[rank] > 0)
879                {
880                        unpackIntersection(elements, recvBuffer2[rank]);
881                        delete [] recvBuffer2[rank];
882                }
883        }
884        delete [] sendBuffer2;
885        delete [] recvBuffer2;
886        delete [] sendBuffer;
887        delete [] recvBuffer;
888
889        delete [] nbSendNode;
890        delete [] nbRecvNode;
891        delete [] sentMessageSize;
892        delete [] recvMessageSize;
893}
894
895Mapper::~Mapper()
896{
897        delete [] remapMatrix;
898        delete [] srcAddress;
899        delete [] srcRank;
900        delete [] dstAddress;
901        if (neighbourElements) delete [] neighbourElements;
902}
903
904}
Note: See TracBrowser for help on using the repository browser.