source: XIOS/dev/branch_openmp/extern/remap/src/parallel_tree.cpp @ 1460

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

branch_openmp merged with XIOS_DEV_CMIP6@1459

File size: 16.5 KB
Line 
1#include <cassert>
2#include "node.hpp"
3#include "timerRemap.hpp"
4#include "circle.hpp"
5#include "meshutil.hpp"
6#include "polyg.hpp"
7#include "gridRemap.hpp"
8#include "intersect.hpp"
9#include "errhandle.hpp"
10#include "mpi_routing.hpp"
11#include "misc.hpp"
12
13#include "parallel_tree.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
24static const int assignLevel = 2;
25
26// only the circle is packed, rest of node will be initialized on unpacking
27static void packNode(Node& node, char *buffer, int& index)
28{
29        if (buffer == NULL) // compute size only
30                index += 4 * sizeof(double);
31        else
32        {
33                *(Coord *)(&buffer[index]) = node.centre;
34                index += sizeof(Coord);
35
36                *(double *)(&buffer[index]) = node.radius;
37                index += sizeof(double);
38        }
39}
40
41static void unpackNode(Node& node, char *buffer, int& index)
42{
43        Coord centre;
44        double r;
45
46        centre = *(Coord *)(&buffer[index]);
47        index += sizeof(Coord);
48
49        r = *(double *)(&buffer[index]);
50        index += sizeof(double);
51
52        node.centre = centre;
53        node.radius = r;
54}
55
56
57static void packElement(Elt *ptElement, char *buffer, int& index)
58{
59        if (buffer == NULL) index += packedPolygonSize(*ptElement);
60        else packPolygon(*ptElement, buffer, index);
61}
62
63static void unpackElement(Elt *ptElement, char *buffer, int& index)
64{
65        unpackPolygon(*ptElement, buffer, index);
66}
67
68static void packVector(const vector<int>& vec, char *buf, int& pos)
69{
70        if (buf == NULL)
71                pos += sizeof(int) + vec.size()*sizeof(int);
72        else
73        {
74                *((int *) &(buf[pos])) = vec.size();
75                pos += sizeof(int);
76                for (int i = 0; i < vec.size(); i++)
77                {
78                        *((int *) &(buf[pos])) = vec[i];
79                        pos += sizeof(int);
80                }
81        }
82}
83
84static void unpackVector(vector<int>& vec, const char *buf, int& pos)
85{
86        vec.resize(*((int *) &(buf[pos])));
87        pos += sizeof(int);
88        for (int i = 0; i < vec.size(); i++)
89        {
90                vec[i] = *((int *) &(buf[pos]));
91                pos += sizeof(int);
92        }
93}
94
95static void assignRoute(CSampleTree& tree, const CCascadeLevel& cl)  // newroot <- root
96{
97        vector<int> routeRank(cl.group_size);
98        for (int i = 0; i < cl.group_size; i++)
99                routeRank[i] = i; //(cl.group_size + cl.polour() - i) % cl.group_size;
100        std::vector<int>::iterator rank = routeRank.begin();
101        tree.root->assignRoute(rank, assignLevel);
102}
103
104static void transferRoutedNodes(CMPIRouting& MPIRoute, /*const*/ vector<Node>& nodeSend, const vector<int>& route, vector<Node>& nodeRecv)
105{
106        /* `route` knows what goes where */
107        MPIRoute.init(route);
108        int nRecv = MPIRoute.getTotalSourceElement();
109        nodeRecv.resize(nRecv);
110        MPIRoute.transferToTarget(&nodeSend[0], &nodeRecv[0], packNode, unpackNode);
111}
112
113static void transferRoutedIntersections(CMPIRouting& MPIRoute, const vector<Node>& nodeSend, const vector<vector<int> >& route, vector<Node>& nodeRecv)
114{
115        // `route` knows what goes where
116        MPIRoute.init(route);
117        int nRecv = MPIRoute.getTotalSourceElement();
118        nodeRecv.resize(nRecv);
119        MPIRoute.transferToTarget((Node * /*mpi wants non-const*/) &nodeSend[0], &nodeRecv[0], packNode, unpackNode);
120//cout << MPIRoute.mpiRank << " ROUTE " << nRecv << ": " << nodeSend.size() << " " << nodeRecv.size() << "    ";
121}
122
123//CParallelTree::CParallelTree(MPI_Comm comm) : communicator(comm), cascade(MIN_NODE_SZ*MIN_NODE_SZ*2, comm)
124CParallelTree::CParallelTree(MPI_Comm comm) : cascade(MAX_NODE_SZ*MAX_NODE_SZ*2, comm)
125{
126  MPI_Comm_dup(comm, &communicator);
127  treeCascade.reserve(cascade.num_levels);
128  for (int lev = 0; lev < cascade.num_levels; lev++)
129    treeCascade.push_back(CSampleTree(cascade.level[lev].group_size, assignLevel));
130}
131
132void CParallelTree::buildSampleTreeCascade(vector<Node>& sampleNodes /*route field will be modified*/, int level)
133{
134        buildSampleTree(treeCascade[level], sampleNodes, cascade.level[level]);
135        assignRoute(treeCascade[level] /*out*/, cascade.level[level] /*in*/);
136
137        if (level+1 < cascade.num_levels)
138        {
139                vector<int> route(sampleNodes.size());
140                treeCascade[level].routeNodes(route, sampleNodes, assignLevel);
141
142                vector<Node> routedNodes;
143                CMPIRouting mpiRoute(cascade.level[level].pg_comm);
144                transferRoutedNodes(mpiRoute, sampleNodes, route, routedNodes);
145                buildSampleTreeCascade(routedNodes, level+1);
146        }
147}
148
149void buildSampleTree(CSampleTree& tree, const vector<Node>& node, const CCascadeLevel& comm)
150/*
151        In the beginning all the sample elements are distributed
152        -> communicate to make available at each rank
153           so that each rank can build the same sample tree
154*/
155{
156        int n = node.size(); // number of samples intially on this rank (local)
157
158        int blocSize = comm.group_size * ipow(MAX_NODE_SZ, assignLevel);
159
160        int nrecv; // global number of samples  THIS WILL BE THE NUMBER OF LEAFS IN THE SAMPLE TREE
161        MPI_Allreduce(&n, &nrecv, 1, MPI_INT, MPI_SUM, comm.comm); // => size of sample tree does not depend on keepNodes!
162
163        double ratio = blocSize / (1.0 * nrecv);
164        int nsend = ratio * n + 1; // nsend = n_local_samples / n_global_samples * blocksize + 1 = blocksize/comm.size
165        if (nsend > n) nsend = n;
166
167        int *counts = new int[comm.size];
168        MPI_Allgather(&nsend, 1, MPI_INT, counts, 1, MPI_INT, comm.comm);
169
170        nrecv = 0;
171        int *displs = new int[comm.size];
172        for (int i = 0; i < comm.size; i++)
173        {
174                displs[i] = 4 * nrecv;
175                nrecv += counts[i];
176                counts[i] = 4 * counts[i];
177        }
178
179        /* pack circles around sample elements */
180        double *sendBuffer = new double[nsend*4];
181        int index = 0;
182        vector<int> randomArray(n);
183        randomizeArray(randomArray);
184        for (int i = 0; i < nsend; i++)
185        {
186                const Node& no = node[randomArray[i]];
187                *((Coord *) (sendBuffer + index)) = no.centre;
188                index += sizeof(Coord)/sizeof(*sendBuffer);
189                sendBuffer[index++] = no.radius;
190        }
191
192        /* each process needs the sample elements from all processes */
193        double *recvBuffer = new double[nrecv*4];
194        MPI_Allgatherv(sendBuffer, 4 * nsend, MPI_DOUBLE, recvBuffer, counts, displs, MPI_DOUBLE, comm.comm);
195        delete[] sendBuffer;
196        delete[] counts;
197        delete[] displs;
198
199        /* unpack */
200/*
201        randomArray.resize(nrecv);
202        randomizeArray(randomArray);
203        tree.leafs.resize(nrecv);
204        index = 0;
205
206
207  for (int i = 0; i < nrecv; i++)
208        {
209                Coord x = *(Coord *)(&recvBuffer[index]);
210                index += sizeof(Coord)/sizeof(*recvBuffer);
211                double radius = recvBuffer[index++];
212                tree.leafs[randomArray[i]].centre = x;
213                tree.leafs[randomArray[i]].radius = radius;
214
215        }
216*/
217
218  randomArray.resize(blocSize);
219        randomizeArray(randomArray);
220        tree.leafs.resize(blocSize);
221        index = 0;
222 
223  size_t s=(sizeof(Coord)/sizeof(*recvBuffer)+1)*nrecv ;
224 
225  for (int i = 0; i < blocSize; i++)
226        {
227                Coord x = *(Coord *)(&recvBuffer[index%s]);
228                index += sizeof(Coord)/sizeof(*recvBuffer);
229                double radius = recvBuffer[index%s];
230    index++ ;
231                tree.leafs[randomArray[i]].centre = x;
232                tree.leafs[randomArray[i]].radius = radius;
233
234        }
235
236
237        delete [] recvBuffer;
238
239        CTimer::get("buildSampleTree(local)").resume();
240        tree.build(tree.leafs);
241//      cout << "SampleTree build : assign Level " << assignLevel << " nb Nodes : " << tree.levelSize[assignLevel] << endl;
242        CTimer::get("buildSampleTree(local)").suspend();
243        CTimer::get("buildSampleTree(local)").print();
244
245        /* End gracefully if sample tree could not be constructed with desired number of nodes on assignLevel */
246        int allok, ok = (tree.levelSize[assignLevel] == comm.group_size);
247        if (!ok)
248  {
249    cerr << comm.rank << ": PROBLEM: (node assign)" << tree.levelSize[assignLevel] << " != " << comm.group_size << " (keepNodes)" 
250         << "   node size : "<<node.size()<<"   bloc size : "<<blocSize<<"  total number of leaf : "<<tree.leafs.size()<<endl ;
251/*
252        MPI_Allreduce(&ok, &allok, 1, MPI_INT, MPI_PROD, communicator);
253        if (!allok) {
254                MPI_Finalize();
255                exit(1);
256        }
257*/
258    MPI_Abort(MPI_COMM_WORLD,-1) ;
259  }
260/*
261  cout<<"*******************************************"<<endl ;
262  cout<<"******* Sample Tree output        *********"<<endl ;
263  cout<<"*******************************************"<<endl ;
264  tree.output(cout,1) ;
265
266  cout<<"*******************************************"<<endl ;
267*/
268
269  assert(tree.root->incluCheck() == 0);
270}
271
272
273void CParallelTree::buildLocalTree(const vector<Node>& node, const vector<int>& route)
274{
275        CMPIRouting MPIRoute(communicator);
276        MPI_Barrier(communicator);
277        CTimer::get("buildLocalTree(initRoute)").resume();
278        MPIRoute.init(route);
279        CTimer::get("buildLocalTree(initRoute)").suspend();
280        CTimer::get("buildLocalTree(initRoute)").print();
281
282        nbLocalElements = MPIRoute.getTotalSourceElement();
283        localElements = new Elt[nbLocalElements];
284
285        vector<Elt*> ptElement(node.size());
286        for (int i = 0; i < node.size(); i++)
287                ptElement[i] = (Elt *) (node[i].data);
288
289        vector<Elt*> ptLocalElement(nbLocalElements);
290        for (int i = 0; i < nbLocalElements; i++)
291                ptLocalElement[i] = &localElements[i];
292
293        CTimer::get("buildLocalTree(transfer)").resume();
294        MPIRoute.transferToTarget(&ptElement[0], &ptLocalElement[0], packElement, unpackElement);
295        CTimer::get("buildLocalTree(transfer)").suspend();
296        CTimer::get("buildLocalTree(transfer)").print();
297
298        CTimer::get("buildLocalTree(local)").resume();
299
300        int mpiRank;
301        MPI_Comm_rank(communicator, &mpiRank);
302        localTree.leafs.reserve(nbLocalElements);
303        for (int i = 0; i < nbLocalElements; i++)
304        {
305                Elt& elt = localElements[i];
306                elt.id.ind = i;
307                elt.id.rank = mpiRank;
308                localTree.leafs.push_back(Node(elt.x, cptRadius(elt), &localElements[i]));
309        }
310        localTree.build(localTree.leafs);
311
312        cptAllEltsGeom(localElements, nbLocalElements, srcGrid.pole);
313        CTimer::get("buildLocalTree(local)").suspend();
314        CTimer::get("buildLocalTree(local)").print();
315}
316
317void CParallelTree::build(vector<Node>& node, vector<Node>& node2)
318{
319
320        int assignLevel = 2;
321        int nbSampleNodes = 2*ipow(MAX_NODE_SZ + 1, assignLevel);
322
323
324  long int nb1, nb2, nb, nbTot ;
325  nb1=node.size() ; nb2=node2.size() ;
326  nb=nb1+nb2 ;
327  MPI_Allreduce(&nb, &nbTot, 1, MPI_LONG, MPI_SUM, communicator) ;
328  int commSize ;
329  MPI_Comm_size(communicator,&commSize) ;
330 
331        // make multiple of two
332        nbSampleNodes /= 2;
333        nbSampleNodes *= 2;
334//  assert( nbTot > nbSampleNodes*commSize) ;
335   
336  int nbSampleNodes1 = nbSampleNodes * (nb1*commSize)/(1.*nbTot) ;
337  int nbSampleNodes2 = nbSampleNodes * (nb2*commSize)/(1.*nbTot) ;
338 
339
340//      assert(node.size() > nbSampleNodes);
341//      assert(node2.size() > nbSampleNodes);
342//      assert(node.size() + node2.size() > nbSampleNodes);
343        vector<Node> sampleNodes; sampleNodes.reserve(nbSampleNodes1+nbSampleNodes2);
344
345        vector<int> randomArray1(node.size());
346        randomizeArray(randomArray1);
347        vector<int> randomArray2(node2.size());
348        randomizeArray(randomArray2);
349
350/*     
351        int s1,s2 ;
352       
353        if (node.size()< nbSampleNodes/2)
354        {
355          s1 = node.size() ;
356          s2 = nbSampleNodes-s1 ;
357        }
358        else if (node2.size()< nbSampleNodes/2)
359        {
360          s2 = node.size() ;
361          s1 = nbSampleNodes-s2 ;
362        }
363        else
364        {
365          s1=nbSampleNodes/2 ;
366          s2=nbSampleNodes/2 ;
367        }
368*/
369        for (int i = 0; i <nbSampleNodes1; i++) sampleNodes.push_back(Node(node[randomArray1[i%nb1]].centre,  node[randomArray1[i%nb1]].radius, NULL));
370        for (int i = 0; i <nbSampleNodes2; i++) sampleNodes.push_back(Node(node2[randomArray2[i%nb2]].centre, node2[randomArray2[i%nb2]].radius, NULL));
371
372/*         
373        for (int i = 0; i < nbSampleNodes/2; i++)
374        {
375          sampleNodes.push_back(Node(node[randomArray1[i]].centre,  node[randomArray1[i]].radius, NULL));
376          sampleNodes.push_back(Node(node2[randomArray2[i]].centre, node2[randomArray2[i]].radius, NULL));
377        }
378*/
379        CTimer::get("buildParallelSampleTree").resume();
380        //sampleTree.buildParallelSampleTree(sampleNodes, cascade);
381        buildSampleTreeCascade(sampleNodes);
382        CTimer::get("buildParallelSampleTree").suspend();
383        CTimer::get("buildParallelSampleTree").print();
384
385        //route source mesh
386        CTimer::get("parallelRouteNode").resume();
387        vector<int> route(node.size());
388        cout<<"node.size = "<<node.size()<<endl;
389        routeNodes(route /*out*/, node);
390        CTimer::get("parallelRouteNode").suspend();
391        CTimer::get("parallelRouteNode").print();
392
393        CTimer::get("buildLocalTree").resume();
394        buildLocalTree(node, route);
395        CTimer::get("buildLocalTree").suspend();
396        CTimer::get("buildLocalTree").print();
397
398        CTimer::get("buildRouteTree").resume();
399        /* update circles of tree cascade so it can be used for routing */
400        updateCirclesForRouting(localTree.root->centre, localTree.root->radius);
401        CTimer::get("buildRouteTree").suspend();
402        CTimer::get("buildRouteTree").print();
403}
404
405void CParallelTree::routeNodes(vector<int>& route, vector<Node>& nodes /*route field used*/, int level)
406{
407        treeCascade[level].routeNodes(route /*assign*/, nodes, assignLevel);
408
409        if (level+1 < cascade.num_levels)
410        {
411                vector<Node> routedNodes;
412                CMPIRouting MPIRoute(cascade.level[level].pg_comm);
413                transferRoutedNodes(MPIRoute, nodes, route /*use*/, routedNodes);
414                vector<int> globalRank(routedNodes.size());
415                routeNodes(globalRank, routedNodes, level + 1);
416                MPIRoute.transferFromSource(&route[0] /*override*/, &globalRank[0]);
417        }
418        else
419        {
420                CMPIRouting MPIRoute(cascade.level[level].comm); // or use pg_comm, on last cascade level identical
421                MPIRoute.init(route);
422                int nbRecvNode = MPIRoute.getTotalSourceElement();
423                vector<int> globalRank(nbRecvNode);
424                for (int i = 0; i < globalRank.size(); i++)
425                        globalRank[i] = cascade.level[0].rank;
426                MPIRoute.transferFromSource(&route[0] /*override*/, &globalRank[0]);
427        }
428}
429
430/* assume `to` to be empty vector at entry */
431void linearize(const vector<vector<int> >& from, vector<int>& to)
432{
433        int cnt = 0;
434        for (int i = 0; i < from.size(); ++i)
435                cnt += from[i].size();
436        to.resize(cnt);
437        vector<int>::iterator dest = to.begin();
438        for (int i = 0; i < from.size(); ++i)
439                dest = copy(from[i].begin(), from[i].end(), dest);
440}
441
442/* at entry `to` must already have it's final shape and only values are overriden */
443void delinearize(const vector<int>& from, vector<vector<int> >& to)
444{
445        vector<int>::const_iterator end, src = from.begin();
446        for (int i = 0; i < to.size(); ++i, src=end)
447                copy(src, end = src + to[i].size(), to[i].begin());
448}
449
450void CParallelTree::routeIntersections(vector<vector<int> >& routes, vector<Node>& nodes, int level)
451{
452        treeCascade[level].routeIntersections(routes /*assign*/, nodes);
453
454        if (level+1 < cascade.num_levels)
455        {
456                vector<Node> routedNodes;
457                CMPIRouting MPIRoute(cascade.level[level].pg_comm);
458
459                vector<int> flattenedRoutes1;
460                linearize(routes, flattenedRoutes1);
461                vector<Node> double_nodes(flattenedRoutes1.size());
462                int j = 0;
463                for (int i = 0; i < routes.size(); ++i)
464                        for (int k = 0; k < routes[i].size(); ++k, ++j)
465                                double_nodes[j] = nodes[i];
466                transferRoutedNodes(MPIRoute, double_nodes, flattenedRoutes1 /*use*/, routedNodes);
467                vector<vector<int> > globalRanks(routedNodes.size());
468                routeIntersections(globalRanks /*out*/, routedNodes /*in*/, level + 1);
469                vector<vector<int> > flattenedRoutes(flattenedRoutes1.size());
470                // transferFromSource expects sizes (nbSendNode=flattenedRoutes, nbRecvNode=routedNodes.size())
471                MPIRoute.transferFromSource(&flattenedRoutes[0], &globalRanks[0], packVector, unpackVector);
472                for (int i = 0, j = 0; i < routes.size(); ++i)
473                {
474                        int old_size = routes[i].size();
475                        routes[i].resize(0);
476                        for (int k = 0; k < old_size; ++k, ++j)
477                                for (int l = 0; l < flattenedRoutes[j].size(); ++l)
478                                        routes[i].push_back(flattenedRoutes[j][l]);
479                }
480                assert(j == flattenedRoutes1.size());
481
482        }
483        else
484        {
485                CMPIRouting MPIRoute(cascade.level[level].comm);
486                MPIRoute.init(routes);
487                int nbRecvNode = MPIRoute.getTotalSourceElement();
488                vector<int> globalRanks(nbRecvNode, cascade.level[0].rank);
489                vector<int> flattened_routes;
490                linearize(routes, flattened_routes);
491                MPIRoute.transferFromSource(&flattened_routes[0], &globalRanks[0]);
492                delinearize(flattened_routes, routes);
493        }
494        if (level!=level)
495        {
496                for (int i = 0; i < routes.size(); ++i)
497                        for (int k = 0; k < routes[i].size(); ++k)
498                                if (routes[i][k] == cascade.level[0].rank) routes[i].erase(routes[i].begin() + (k--));
499        }
500}
501
502void CParallelTree::updateCirclesForRouting(Coord rootCentre, double rootRadius, int level)
503{
504        if (level + 1 < cascade.num_levels) // children in the cascade have to update first
505        {
506                updateCirclesForRouting(rootCentre, rootRadius, level + 1);
507                rootCentre = treeCascade[level+1].root->centre;
508                rootRadius = treeCascade[level+1].root->radius;
509        }
510
511        // gather circles on this level of the cascade
512        int pg_size;
513        MPI_Comm_size(cascade.level[level].pg_comm, &pg_size);
514        vector<Coord> allRootCentres(pg_size);
515        vector<double> allRootRadia(pg_size);
516        MPI_Allgather(&rootCentre, 3, MPI_DOUBLE, &allRootCentres[0], 3, MPI_DOUBLE, cascade.level[level].pg_comm);
517        MPI_Allgather(&rootRadius, 1, MPI_DOUBLE, &allRootRadia[0],   1, MPI_DOUBLE, cascade.level[level].pg_comm);
518
519        // now allRootsRadia and allRootCentres must be inserted into second levels of us and propagated to root
520        treeCascade[level].root->assignCircleAndPropagateUp(&allRootCentres[0], &allRootRadia[0], assignLevel);
521}
522
523CParallelTree::~CParallelTree()
524{
525        delete [] localElements;
526}
527
528}
Note: See TracBrowser for help on using the repository browser.