source: XIOS/dev/branch_yushan_merged/extern/remap/src/mapper.cpp @ 1149

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

bug fixed in mpi_comm_split. Key needs to be specifify.

File size: 32.9 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, bool renormalize, bool quantity)
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, renormalize, quantity);
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, bool renormalize, bool quantity)
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    double **recvArea = 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            recvArea[rank]    = new double[nbSendElement[rank]];
209            if (order == 2)
210            {
211                recvNeighIds[rank] = new GloId[nbSendElement[rank]*(NMAX+1)];
212                recvGrad[rank]    = new Coord[nbSendElement[rank]*(NMAX+1)];
213            }
214            else
215                recvNeighIds[rank] = new GloId[nbSendElement[rank]];
216
217            last = -1;
218            index = -1;
219            for (multimap<int, Polyg *>::iterator it = elementList[rank].begin(); it != elementList[rank].end(); ++it)
220            {
221                if (last != it->first)
222                    index++;
223                sendElement[rank][index] = it->first;
224                last = it->first;
225            }
226        }
227    }
228
229    /* communicate sizes of source elements to be sent (index lists and later values and gradients) */
230    int *nbRecvElement = new int[mpiSize];
231    MPI_Alltoall(nbSendElement, 1, MPI_INT, nbRecvElement, 1, MPI_INT, communicator);
232
233    /* communicate indices of source elements on other ranks whoes value and gradient we need (since intersection) */
234    int nbSendRequest = 0;
235    int nbRecvRequest = 0;
236    int **recvElement = new int*[mpiSize];
237    double **sendValue = new double*[mpiSize];
238    double **sendArea = new double*[mpiSize];
239    Coord **sendGrad = new Coord*[mpiSize];
240    GloId **sendNeighIds = new GloId*[mpiSize];
241    MPI_Request *sendRequest = new MPI_Request[4*mpiSize];
242    MPI_Request *recvRequest = new MPI_Request[4*mpiSize];
243    for (int rank = 0; rank < mpiSize; rank++)
244    {
245        if (nbSendElement[rank] > 0)
246        {
247            MPI_Issend(sendElement[rank], nbSendElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
248            nbSendRequest++;
249        }
250
251        if (nbRecvElement[rank] > 0)
252        {
253            recvElement[rank] = new int[nbRecvElement[rank]];
254            sendValue[rank]   = new double[nbRecvElement[rank]];
255            sendArea[rank]   = new double[nbRecvElement[rank]];
256            if (order == 2)
257            {
258                sendNeighIds[rank] = new GloId[nbRecvElement[rank]*(NMAX+1)];
259                sendGrad[rank]    = new Coord[nbRecvElement[rank]*(NMAX+1)];
260            }
261            else
262            {
263                sendNeighIds[rank] = new GloId[nbRecvElement[rank]];
264            }
265            MPI_Irecv(recvElement[rank], nbRecvElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
266            nbRecvRequest++;
267        }
268    }
269
270    MPI_Status *status = new MPI_Status[4*mpiSize];
271   
272    MPI_Waitall(nbRecvRequest, recvRequest, status);
273    MPI_Waitall(nbSendRequest, sendRequest, status);
274   
275
276    /* for all indices that have been received from requesting ranks: pack values and gradients, then send */
277    nbSendRequest = 0;
278    nbRecvRequest = 0;
279    for (int rank = 0; rank < mpiSize; rank++)
280    {
281        if (nbRecvElement[rank] > 0)
282        {
283            int jj = 0; // jj == j if no weight writing
284            for (int j = 0; j < nbRecvElement[rank]; j++)
285            {
286                sendValue[rank][j] = sstree.localElements[recvElement[rank][j]].val;
287                sendArea[rank][j] = sstree.localElements[recvElement[rank][j]].area;
288                if (order == 2)
289                {
290                    sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].grad;
291                    sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].src_id;
292                    jj++;
293                    for (int i = 0; i < NMAX; i++)
294                    {
295                        sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].gradNeigh[i];
296                        sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].neighId[i];
297                        jj++;
298                    }
299                }
300                else
301                    sendNeighIds[rank][j] = sstree.localElements[recvElement[rank][j]].src_id;
302            }
303            MPI_Issend(sendValue[rank], nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]);
304            nbSendRequest++;
305            MPI_Issend(sendArea[rank], nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]);
306            nbSendRequest++;
307            if (order == 2)
308            {
309                MPI_Issend(sendGrad[rank], 3*nbRecvElement[rank]*(NMAX+1),
310                        MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]);
311                nbSendRequest++;
312                MPI_Issend(sendNeighIds[rank], 4*nbRecvElement[rank]*(NMAX+1), MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
313                //ym  --> attention taille GloId
314                nbSendRequest++;
315            }
316            else
317            {
318                MPI_Issend(sendNeighIds[rank], 4*nbRecvElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]);
319                //ym  --> attention taille GloId
320                nbSendRequest++;               
321            }
322        }
323        if (nbSendElement[rank] > 0)
324        {
325            MPI_Irecv(recvValue[rank],  nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]);
326            nbRecvRequest++;
327            MPI_Irecv(recvArea[rank],  nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]);
328            nbRecvRequest++;
329            if (order == 2)
330            {
331                MPI_Irecv(recvGrad[rank], 3*nbSendElement[rank]*(NMAX+1),
332                        MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]);
333                nbRecvRequest++;
334                MPI_Irecv(recvNeighIds[rank], 4*nbSendElement[rank]*(NMAX+1), MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
335                //ym  --> attention taille GloId
336                nbRecvRequest++;
337            }
338            else
339            {
340                MPI_Irecv(recvNeighIds[rank], 4*nbSendElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]);
341                //ym  --> attention taille GloId
342                nbRecvRequest++;
343            }
344        }
345    }
346   
347    MPI_Waitall(nbRecvRequest, recvRequest, status);
348    MPI_Waitall(nbSendRequest, sendRequest, status); 
349   
350   
351
352    /* now that all values and gradients are available use them to computed interpolated values on target
353       and also to compute weights */
354    int i = 0;
355    for (int j = 0; j < nbElements; j++)
356    {
357        Elt& e = elements[j];
358
359        /* since for the 2nd order case source grid elements can contribute to a destination grid element over several "paths"
360           (step1: gradient is computed using neighbours on same grid, step2: intersection uses several elements on other grid)
361           accumulate them so that there is only one final weight between two elements */
362        map<GloId,double> wgt_map;
363
364        /* for destination element `e` loop over all intersetions/the corresponding source elements */
365        for (list<Polyg *>::iterator it = e.is.begin(); it != e.is.end(); it++)
366        {
367            /* it is the intersection element, so it->x and it->area are barycentre and area of intersection element (super mesh)
368               but it->id is id of the source element that it intersects */
369            int n1 = (*it)->id.ind;
370            int rank = (*it)->id.rank;
371            double fk = recvValue[rank][n1];
372            double srcArea = recvArea[rank][n1];
373            double w = (*it)->area;
374            if (quantity) w/=srcArea ;
375
376            /* first order: src value times weight (weight = supermesh area), later divide by target area */
377            int kk = (order == 2) ? n1 * (NMAX + 1) : n1;
378            GloId neighID = recvNeighIds[rank][kk];
379            wgt_map[neighID] += w;
380
381            if (order == 2)
382            {
383                for (int k = 0; k < NMAX+1; k++)
384                {
385                    int kk = n1 * (NMAX + 1) + k;
386                    GloId neighID = recvNeighIds[rank][kk];
387                    if (neighID.ind != -1)  wgt_map[neighID] += w * scalarprod(recvGrad[rank][kk], (*it)->x);
388                }
389
390            }
391        }
392
393        double renorm=0;
394        if (renormalize) 
395            for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++) renorm+=it->second / e.area;
396        else renorm=1. ;
397
398        for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++)
399        {
400            if (quantity)  this->remapMatrix[i] = (it->second ) / renorm;
401            else this->remapMatrix[i] = (it->second / e.area) / renorm;
402            this->srcAddress[i] = it->first.ind;
403            this->srcRank[i] = it->first.rank;
404            this->dstAddress[i] = j;
405            this->sourceWeightId[i]= it->first.globalId ;
406            this->targetWeightId[i]= targetGlobalId[j] ;
407            i++;
408        }
409    }
410   
411    /* free all memory allocated in this function */
412    /*for (int rank = 0; rank < mpiSize; rank++)
413    {
414        if (nbSendElement[rank] > 0)
415        {
416            delete[] sendElement[rank];
417            delete[] recvValue[rank];
418            delete[] recvArea[rank];
419            if (order == 2)
420            {
421                delete[] recvGrad[rank];
422            }
423            delete[] recvNeighIds[rank];
424        }
425        if (nbRecvElement[rank] > 0)
426        {
427            delete[] recvElement[rank];
428            delete[] sendValue[rank];
429            delete[] sendArea[rank];
430            if (order == 2)
431                delete[] sendGrad[rank];
432            delete[] sendNeighIds[rank];
433        }
434    }
435    delete[] status;
436    delete[] sendRequest;
437    delete[] recvRequest;
438    delete[] elementList;
439    delete[] nbSendElement;
440    delete[] nbRecvElement;
441    delete[] sendElement;
442    delete[] recvElement;
443    delete[] sendValue;
444    delete[] recvValue;
445    delete[] sendGrad;
446    delete[] recvGrad;
447    delete[] sendNeighIds;
448    delete[] recvNeighIds;
449    */
450    return i;
451}
452
453void Mapper::computeGrads()
454{
455    /* array of pointers to collect local elements and elements received from other cpu */
456    vector<Elt*> globalElements(sstree.nbLocalElements + nbNeighbourElements);
457    int index = 0;
458    for (int i = 0; i < sstree.nbLocalElements; i++, index++)
459        globalElements[index] = &(sstree.localElements[i]);
460    for (int i = 0; i < nbNeighbourElements; i++, index++)
461        globalElements[index] = &neighbourElements[i];
462
463    update_baryc(sstree.localElements, sstree.nbLocalElements);
464    computeGradients(&globalElements[0], sstree.nbLocalElements);
465}
466
467/** for each element of the source grid, finds all the neighbouring elements that share an edge
468  (filling array neighbourElements). This is used later to compute gradients */
469void Mapper::buildMeshTopology()
470{
471    int mpiSize, mpiRank;
472    MPI_Comm_size(communicator, &mpiSize);
473    MPI_Comm_rank(communicator, &mpiRank);
474
475    vector<Node> *routingList = new vector<Node>[mpiSize];
476    vector<vector<int> > routes(sstree.localTree.leafs.size());
477
478    sstree.routeIntersections(routes, sstree.localTree.leafs);
479
480    for (int i = 0; i < routes.size(); ++i)
481        for (int k = 0; k < routes[i].size(); ++k)
482            routingList[routes[i][k]].push_back(sstree.localTree.leafs[i]);
483    routingList[mpiRank].clear();
484
485
486    CMPIRouting mpiRoute(communicator);
487    mpiRoute.init(routes);
488    int nRecv = mpiRoute.getTotalSourceElement();
489
490    int *nbSendNode = new int[mpiSize];
491    int *nbRecvNode = new int[mpiSize];
492    int *sendMessageSize = new int[mpiSize];
493    int *recvMessageSize = new int[mpiSize];
494
495    for (int rank = 0; rank < mpiSize; rank++)
496    {
497        nbSendNode[rank] = routingList[rank].size();
498        sendMessageSize[rank] = 0;
499        for (size_t j = 0; j < routingList[rank].size(); j++)
500        {
501            Elt *elt = (Elt *) (routingList[rank][j].data);
502            sendMessageSize[rank] += packedPolygonSize(*elt);
503        }
504    }
505
506    MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
507    MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
508
509    char **sendBuffer = new char*[mpiSize];
510    char **recvBuffer = new char*[mpiSize];
511    int *pos = new int[mpiSize];
512
513    for (int rank = 0; rank < mpiSize; rank++)
514    {
515        if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sendMessageSize[rank]];
516        if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]];
517    }
518
519    for (int rank = 0; rank < mpiSize; rank++)
520    {
521        pos[rank] = 0;
522        for (size_t j = 0; j < routingList[rank].size(); j++)
523        {
524            Elt *elt = (Elt *) (routingList[rank][j].data);
525            packPolygon(*elt, sendBuffer[rank], pos[rank]);
526        }
527    }
528    delete [] routingList;
529
530
531    int nbSendRequest = 0;
532    int nbRecvRequest = 0;
533    MPI_Request *sendRequest = new MPI_Request[mpiSize];
534    MPI_Request *recvRequest = new MPI_Request[mpiSize];
535    MPI_Status  *status      = new MPI_Status[mpiSize];
536
537    for (int rank = 0; rank < mpiSize; rank++)
538    {
539        if (nbSendNode[rank] > 0)
540        {
541            MPI_Issend(sendBuffer[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
542            nbSendRequest++;
543        }
544        if (nbRecvNode[rank] > 0)
545        {
546            MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
547            nbRecvRequest++;
548        }
549    }
550   
551    MPI_Waitall(nbRecvRequest, recvRequest, status);
552    MPI_Waitall(nbSendRequest, sendRequest, status);
553
554    for (int rank = 0; rank < mpiSize; rank++)
555        if (nbSendNode[rank] > 0) delete [] sendBuffer[rank];
556    delete [] sendBuffer;
557
558    char **sendBuffer2 = new char*[mpiSize];
559    char **recvBuffer2 = new char*[mpiSize];
560
561    for (int rank = 0; rank < mpiSize; rank++)
562    {
563        nbSendNode[rank] = 0;
564        sendMessageSize[rank] = 0;
565
566        if (nbRecvNode[rank] > 0)
567        {
568            set<NodePtr> neighbourList;
569            pos[rank] = 0;
570            for (int j = 0; j < nbRecvNode[rank]; j++)
571            {
572                Elt elt;
573                unpackPolygon(elt, recvBuffer[rank], pos[rank]);
574                Node node(elt.x, cptRadius(elt), &elt);
575                findNeighbour(sstree.localTree.root, &node, neighbourList);
576            }
577            nbSendNode[rank] = neighbourList.size();
578            for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++)
579            {
580                Elt *elt = (Elt *) ((*it)->data);
581                sendMessageSize[rank] += packedPolygonSize(*elt);
582            }
583
584            sendBuffer2[rank] = new char[sendMessageSize[rank]];
585            pos[rank] = 0;
586
587            for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++)
588            {
589                Elt *elt = (Elt *) ((*it)->data);
590                packPolygon(*elt, sendBuffer2[rank], pos[rank]);
591            }
592        }
593    }
594    for (int rank = 0; rank < mpiSize; rank++)
595        if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank];
596    delete [] recvBuffer;
597
598
599    MPI_Barrier(communicator);
600    MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
601    MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
602
603    for (int rank = 0; rank < mpiSize; rank++)
604        if (nbRecvNode[rank] > 0) recvBuffer2[rank] = new char[recvMessageSize[rank]];
605
606    nbSendRequest = 0;
607    nbRecvRequest = 0;
608
609    for (int rank = 0; rank < mpiSize; rank++)
610    {
611        if (nbSendNode[rank] > 0)
612        {
613            MPI_Issend(sendBuffer2[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
614            nbSendRequest++;
615        }
616        if (nbRecvNode[rank] > 0)
617        {
618            MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
619            nbRecvRequest++;
620        }
621    }
622
623    MPI_Waitall(nbRecvRequest, recvRequest, status);
624    MPI_Waitall(nbSendRequest, sendRequest, status);
625 
626    int nbNeighbourNodes = 0;
627    for (int rank = 0; rank < mpiSize; rank++)
628        nbNeighbourNodes += nbRecvNode[rank];
629
630    neighbourElements = new Elt[nbNeighbourNodes];
631    nbNeighbourElements = nbNeighbourNodes;
632
633    int index = 0;
634    for (int rank = 0; rank < mpiSize; rank++)
635    {
636        pos[rank] = 0;
637        for (int j = 0; j < nbRecvNode[rank]; j++)
638        {
639            unpackPolygon(neighbourElements[index], recvBuffer2[rank], pos[rank]);
640            neighbourElements[index].id.ind = sstree.localTree.leafs.size() + index;
641            index++;
642        }
643    }
644    for (int rank = 0; rank < mpiSize; rank++)
645    {
646        if (nbRecvNode[rank] > 0) delete [] recvBuffer2[rank];
647        if (nbSendNode[rank] > 0) delete [] sendBuffer2[rank];
648    }
649    delete [] recvBuffer2;
650    delete [] sendBuffer2;
651    delete [] sendMessageSize;
652    delete [] recvMessageSize;
653    delete [] nbSendNode;
654    delete [] nbRecvNode;
655    delete [] sendRequest;
656    delete [] recvRequest;
657    delete [] status;
658    delete [] pos;
659
660    /* re-compute on received elements to avoid having to send this information */
661    neighbourNodes.resize(nbNeighbourNodes);
662    setCirclesAndLinks(neighbourElements, neighbourNodes);
663    cptAllEltsGeom(neighbourElements, nbNeighbourNodes, srcGrid.pole);
664
665    /* the local SS tree must include nodes from other cpus if they are potential
666       intersector of a local node */
667    sstree.localTree.insertNodes(neighbourNodes);
668
669    /* for every local element,
670       use the SS-tree to find all elements (including neighbourElements)
671       who are potential neighbours because their circles intersect,
672       then check all canditates for common edges to build up connectivity information
673       */
674    for (int j = 0; j < sstree.localTree.leafs.size(); j++)
675    {
676        Node& node = sstree.localTree.leafs[j];
677
678        /* find all leafs whoes circles that intersect node's circle and save into node->intersectors */
679        node.search(sstree.localTree.root);
680
681        Elt *elt = (Elt *)(node.data);
682
683        for (int i = 0; i < elt->n; i++) elt->neighbour[i] = NOT_FOUND;
684
685        /* for element `elt` loop through all nodes in the SS-tree
686           whoes circles intersect with the circle around `elt` (the SS intersectors)
687           and check if they are neighbours in the sense that the two elements share an edge.
688           If they do, save this information for elt */
689        for (list<NodePtr>::iterator it = (node.intersectors).begin(); it != (node.intersectors).end(); ++it)
690        {
691            Elt *elt2 = (Elt *)((*it)->data);
692            set_neighbour(*elt, *elt2);
693        }
694
695        /*
696           for (int i = 0; i < elt->n; i++)
697           {
698           if (elt->neighbour[i] == NOT_FOUND)
699           error_exit("neighbour not found");
700           }
701           */
702    }
703}
704
705/** @param elements are the target grid elements */
706void Mapper::computeIntersection(Elt *elements, int nbElements)
707{
708    int mpiSize, mpiRank;
709    MPI_Comm_size(communicator, &mpiSize);
710    MPI_Comm_rank(communicator, &mpiRank);
711
712    MPI_Barrier(communicator);
713
714    vector<Node> *routingList = new vector<Node>[mpiSize];
715
716    vector<Node> routeNodes;  routeNodes.reserve(nbElements);
717    for (int j = 0; j < nbElements; j++)
718    {
719        elements[j].id.ind = j;
720        elements[j].id.rank = mpiRank;
721        routeNodes.push_back(Node(elements[j].x, cptRadius(elements[j]), &elements[j]));
722    }
723
724    vector<vector<int> > routes(routeNodes.size());
725    sstree.routeIntersections(routes, routeNodes);
726    for (int i = 0; i < routes.size(); ++i)
727        for (int k = 0; k < routes[i].size(); ++k)
728            routingList[routes[i][k]].push_back(routeNodes[i]);
729
730    if (verbose >= 2)
731    {
732        cout << " --> rank  " << mpiRank << " nbElements " << nbElements << " : ";
733        for (int rank = 0; rank < mpiSize; rank++)
734            cout << routingList[rank].size() << "   ";
735        cout << endl;
736    }
737    MPI_Barrier(communicator);
738
739    int *nbSendNode = new int[mpiSize];
740    int *nbRecvNode = new int[mpiSize];
741    int *sentMessageSize = new int[mpiSize];
742    int *recvMessageSize = new int[mpiSize];
743
744    for (int rank = 0; rank < mpiSize; rank++)
745    {
746        nbSendNode[rank] = routingList[rank].size();
747        sentMessageSize[rank] = 0;
748        for (size_t j = 0; j < routingList[rank].size(); j++)
749        {
750            Elt *elt = (Elt *) (routingList[rank][j].data);
751            sentMessageSize[rank] += packedPolygonSize(*elt);
752        }
753    }
754
755    MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator);
756    MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
757
758    int total = 0;
759
760    for (int rank = 0; rank < mpiSize; rank++)
761    {
762        total = total + nbRecvNode[rank];
763    }
764
765    if (verbose >= 2) cout << "---> rank " << mpiRank << " : compute intersection : total received nodes  " << total << endl;
766    char **sendBuffer = new char*[mpiSize];
767    char **recvBuffer = new char*[mpiSize];
768    int *pos = new int[mpiSize];
769
770    for (int rank = 0; rank < mpiSize; rank++)
771    {
772        if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sentMessageSize[rank]];
773        if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]];
774    }
775
776    for (int rank = 0; rank < mpiSize; rank++)
777    {
778        pos[rank] = 0;
779        for (size_t j = 0; j < routingList[rank].size(); j++)
780        {
781            Elt* elt = (Elt *) (routingList[rank][j].data);
782            packPolygon(*elt, sendBuffer[rank], pos[rank]);
783        }
784    }
785    delete [] routingList;
786
787    int nbSendRequest = 0;
788    int nbRecvRequest = 0;
789    MPI_Request *sendRequest = new MPI_Request[mpiSize];
790    MPI_Request *recvRequest = new MPI_Request[mpiSize];
791    MPI_Status   *status = new MPI_Status[mpiSize];
792
793    for (int rank = 0; rank < mpiSize; rank++)
794    {
795        if (nbSendNode[rank] > 0)
796        {
797            MPI_Issend(sendBuffer[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
798            nbSendRequest++;
799        }
800        if (nbRecvNode[rank] > 0)
801        {
802            MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
803            nbRecvRequest++;
804        }
805    }
806 
807    MPI_Waitall(nbRecvRequest, recvRequest, status);
808    MPI_Waitall(nbSendRequest, sendRequest, status);
809   
810   
811    char **sendBuffer2 = new char*[mpiSize];
812    char **recvBuffer2 = new char*[mpiSize];
813
814    double tic = cputime();
815    for (int rank = 0; rank < mpiSize; rank++)
816    {
817        sentMessageSize[rank] = 0;
818
819        if (nbRecvNode[rank] > 0)
820        {
821            Elt *recvElt = new Elt[nbRecvNode[rank]];
822            pos[rank] = 0;
823            for (int j = 0; j < nbRecvNode[rank]; j++)
824            {
825                unpackPolygon(recvElt[j], recvBuffer[rank], pos[rank]);
826                cptEltGeom(recvElt[j], tgtGrid.pole);
827                Node recvNode(recvElt[j].x, cptRadius(recvElt[j]), &recvElt[j]);
828                recvNode.search(sstree.localTree.root);
829                /* for a node holding an element of the target, loop throught candidates for intersecting source */
830                for (list<NodePtr>::iterator it = (recvNode.intersectors).begin(); it != (recvNode.intersectors).end(); ++it)
831                {
832                    Elt *elt2 = (Elt *) ((*it)->data);
833                    /* recvElt is target, elt2 is source */
834                    //                                  intersect(&recvElt[j], elt2);
835                    intersect_ym(&recvElt[j], elt2);
836                }
837
838                if (recvElt[j].is.size() > 0) sentMessageSize[rank] += packIntersectionSize(recvElt[j]);
839
840                // here recvNode goes out of scope
841            }
842
843            if (sentMessageSize[rank] > 0)
844            {
845                sentMessageSize[rank] += sizeof(int);
846                sendBuffer2[rank] = new char[sentMessageSize[rank]];
847                *((int *) sendBuffer2[rank]) = 0;
848                pos[rank] = sizeof(int);
849                for (int j = 0; j < nbRecvNode[rank]; j++)
850                {
851                    packIntersection(recvElt[j], sendBuffer2[rank], pos[rank]);
852                    //FIXME should be deleted: recvElt[j].delete_intersections(); // intersection areas have been packed to buffer and won't be used any more
853                }
854            }
855            delete [] recvElt;
856
857        }
858    }
859    delete [] pos;
860
861    for (int rank = 0; rank < mpiSize; rank++)
862    {
863        if (nbSendNode[rank] > 0) delete [] sendBuffer[rank];
864        if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank];
865        nbSendNode[rank] = 0;
866    }
867
868    if (verbose >= 2) cout << "Rank " << mpiRank << "  Compute (internal) intersection " << cputime() - tic << " s" << endl;
869    MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator);
870
871    for (int rank = 0; rank < mpiSize; rank++)
872        if (recvMessageSize[rank] > 0)
873            recvBuffer2[rank] = new char[recvMessageSize[rank]];
874
875    nbSendRequest = 0;
876    nbRecvRequest = 0;
877
878    for (int rank = 0; rank < mpiSize; rank++)
879    {
880        if (sentMessageSize[rank] > 0)
881        {
882            MPI_Issend(sendBuffer2[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]);
883            nbSendRequest++;
884        }
885        if (recvMessageSize[rank] > 0)
886        {
887            MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]);
888            nbRecvRequest++;
889        }
890    }
891   
892    MPI_Waitall(nbRecvRequest, recvRequest, status);
893    MPI_Waitall(nbSendRequest, sendRequest, status);
894   
895   
896
897    delete [] sendRequest;
898    delete [] recvRequest;
899    delete [] status;
900    for (int rank = 0; rank < mpiSize; rank++)
901    {
902        if (nbRecvNode[rank] > 0)
903        {
904            if (sentMessageSize[rank] > 0)
905                delete [] sendBuffer2[rank];
906        }
907
908        if (recvMessageSize[rank] > 0)
909        {
910            unpackIntersection(elements, recvBuffer2[rank]);
911            delete [] recvBuffer2[rank];
912        }
913    }
914    delete [] sendBuffer2;
915    delete [] recvBuffer2;
916    delete [] sendBuffer;
917    delete [] recvBuffer;
918
919    delete [] nbSendNode;
920    delete [] nbRecvNode;
921    delete [] sentMessageSize;
922    delete [] recvMessageSize;
923}
924
925Mapper::~Mapper()
926{
927    delete [] remapMatrix;
928    delete [] srcAddress;
929    delete [] srcRank;
930    delete [] dstAddress;
931    if (neighbourElements) delete [] neighbourElements;
932}
933
934}
Note: See TracBrowser for help on using the repository browser.