Ticket #97: parallel_tree.cpp

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