source: XIOS/dev/branch_openmp/extern/remap/src/mapper.cpp @ 1335

Last change on this file since 1335 was 1335, checked in by yushan, 6 years ago

dev_omp

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