source: XIOS/dev/branch_yushan/extern/remap/src/mapper.cpp @ 1085

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

ep_lib namespace specified when netcdf involved

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