source: XIOS/dev/branch_yushan/extern/remap/src/parallel_tree.cpp @ 1072

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

Using threads : modif for xios_initialize

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