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

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

branch merged with trunk @1200

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