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

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

dev_omp

File size: 15.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"
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, comm)
124CParallelTree::CParallelTree(MPI_Comm comm) : communicator(comm), cascade(MAX_NODE_SZ*MAX_NODE_SZ*2, comm)
125{
126        treeCascade.reserve(cascade.num_levels);
127        for (int lev = 0; lev < cascade.num_levels; lev++)
128                treeCascade.push_back(CSampleTree(cascade.level[lev].group_size, assignLevel));
129}
130
131void CParallelTree::buildSampleTreeCascade(vector<Node>& sampleNodes /*route field will be modified*/, int level)
132{
133        buildSampleTree(treeCascade[level], sampleNodes, cascade.level[level]);
134        assignRoute(treeCascade[level] /*out*/, cascade.level[level] /*in*/);
135
136        if (level+1 < cascade.num_levels)
137        {
138                vector<int> route(sampleNodes.size());
139                treeCascade[level].routeNodes(route, sampleNodes, assignLevel);
140
141                vector<Node> routedNodes;
142                CMPIRouting mpiRoute(cascade.level[level].pg_comm);
143                transferRoutedNodes(mpiRoute, sampleNodes, route, routedNodes);
144                buildSampleTreeCascade(routedNodes, level+1);
145        }
146}
147
148void buildSampleTree(CSampleTree& tree, const vector<Node>& node, const CCascadeLevel& comm)
149/*
150        In the beginning all the sample elements are distributed
151        -> communicate to make available at each rank
152           so that each rank can build the same sample tree
153*/
154{
155        int n = node.size(); // number of samples intially on this rank (local)
156
157        int blocSize = comm.group_size * ipow(MAX_NODE_SZ, assignLevel);
158
159        int nrecv; // global number of samples  THIS WILL BE THE NUMBER OF LEAFS IN THE SAMPLE TREE
160        MPI_Allreduce(&n, &nrecv, 1, MPI_INT, MPI_SUM, comm.comm); // => size of sample tree does not depend on keepNodes!
161
162        double ratio = blocSize / (1.0 * nrecv);
163        int nsend = ratio * n + 1; // nsend = n_local_samples / n_global_samples * blocksize + 1 = blocksize/comm.size
164        if (nsend > n) nsend = n;
165
166        int *counts = new int[comm.size];
167        MPI_Allgather(&nsend, 1, MPI_INT, counts, 1, MPI_INT, comm.comm);
168
169        nrecv = 0;
170        int *displs = new int[comm.size];
171        for (int i = 0; i < comm.size; i++)
172        {
173                displs[i] = 4 * nrecv;
174                nrecv += counts[i];
175                counts[i] = 4 * counts[i];
176        }
177
178        /* pack circles around sample elements */
179        double *sendBuffer = new double[nsend*4];
180        int index = 0;
181        vector<int> randomArray(n);
182        randomizeArray(randomArray);
183        for (int i = 0; i < nsend; i++)
184        {
185                const Node& no = node[randomArray[i]];
186                *((Coord *) (sendBuffer + index)) = no.centre;
187                index += sizeof(Coord)/sizeof(*sendBuffer);
188                sendBuffer[index++] = no.radius;
189        }
190
191        /* each process needs the sample elements from all processes */
192        double *recvBuffer = new double[nrecv*4];
193        MPI_Allgatherv(sendBuffer, 4 * nsend, MPI_DOUBLE, recvBuffer, counts, displs, MPI_DOUBLE, comm.comm);
194        delete[] sendBuffer;
195        delete[] counts;
196        delete[] displs;
197
198
199  randomArray.resize(blocSize);
200        randomizeArray(randomArray);
201        tree.leafs.resize(blocSize);
202        index = 0;
203 
204  size_t s=(sizeof(Coord)/sizeof(*recvBuffer)+1)*nrecv ;
205 
206  for (int i = 0; i < blocSize; i++)
207        {
208                Coord x = *(Coord *)(&recvBuffer[index%s]);
209                index += sizeof(Coord)/sizeof(*recvBuffer);
210                double radius = recvBuffer[index%s];
211    index++ ;
212                tree.leafs[randomArray[i]].centre = x;
213                tree.leafs[randomArray[i]].radius = radius;
214
215        }
216
217
218        delete [] recvBuffer;
219
220        CTimer::get("buildSampleTree(local)").resume();
221        tree.build(tree.leafs);
222//      cout << "SampleTree build : assign Level " << assignLevel << " nb Nodes : " << tree.levelSize[assignLevel] << endl;
223        CTimer::get("buildSampleTree(local)").suspend();
224        CTimer::get("buildSampleTree(local)").print();
225
226        /* End gracefully if sample tree could not be constructed with desired number of nodes on assignLevel */
227        int allok, ok = (tree.levelSize[assignLevel] == comm.group_size);
228        if (!ok)
229  {
230    cerr << comm.rank << ": PROBLEM: (node assign)" << tree.levelSize[assignLevel] << " != " << comm.group_size << " (keepNodes)" 
231         << "   node size : "<<node.size()<<"   bloc size : "<<blocSize<<"  total number of leaf : "<<tree.leafs.size()<<endl ;
232
233    MPI_Abort(MPI_COMM_WORLD,-1) ;
234  }
235/*
236  cout<<"*******************************************"<<endl ;
237  cout<<"******* Sample Tree output        *********"<<endl ;
238  cout<<"*******************************************"<<endl ;
239  tree.output(cout,1) ;
240
241  cout<<"*******************************************"<<endl ;
242*/
243
244  assert(tree.root->incluCheck() == 0);
245}
246
247
248void CParallelTree::buildLocalTree(const vector<Node>& node, const vector<int>& route)
249{
250        CMPIRouting MPIRoute(communicator);
251        MPI_Barrier(communicator);
252        CTimer::get("buildLocalTree(initRoute)").resume();
253        MPIRoute.init(route);
254        CTimer::get("buildLocalTree(initRoute)").suspend();
255        CTimer::get("buildLocalTree(initRoute)").print();
256
257        nbLocalElements = MPIRoute.getTotalSourceElement();
258        localElements = new Elt[nbLocalElements];
259
260        vector<Elt*> ptElement(node.size());
261        for (int i = 0; i < node.size(); i++)
262                ptElement[i] = (Elt *) (node[i].data);
263
264        vector<Elt*> ptLocalElement(nbLocalElements);
265        for (int i = 0; i < nbLocalElements; i++)
266                ptLocalElement[i] = &localElements[i];
267
268        CTimer::get("buildLocalTree(transfer)").resume();
269        MPIRoute.transferToTarget(&ptElement[0], &ptLocalElement[0], packElement, unpackElement);
270        CTimer::get("buildLocalTree(transfer)").suspend();
271        CTimer::get("buildLocalTree(transfer)").print();
272
273        CTimer::get("buildLocalTree(local)").resume();
274
275        int mpiRank;
276        MPI_Comm_rank(communicator, &mpiRank);
277        localTree.leafs.reserve(nbLocalElements);
278        for (int i = 0; i < nbLocalElements; i++)
279        {
280                Elt& elt = localElements[i];
281                elt.id.ind = i;
282                elt.id.rank = mpiRank;
283                localTree.leafs.push_back(Node(elt.x, cptRadius(elt), &localElements[i]));
284        }
285        localTree.build(localTree.leafs);
286
287        cptAllEltsGeom(localElements, nbLocalElements, srcGrid.pole);
288        CTimer::get("buildLocalTree(local)").suspend();
289        CTimer::get("buildLocalTree(local)").print();
290}
291
292void CParallelTree::build(vector<Node>& node, vector<Node>& node2)
293{
294
295        int assignLevel = 2;
296        int nbSampleNodes = 2*ipow(MAX_NODE_SZ + 1, assignLevel);
297
298
299  long int nb1, nb2, nb, nbTot ;
300  nb1=node.size() ; nb2=node2.size() ;
301  nb=nb1+nb2 ;
302  MPI_Allreduce(&nb, &nbTot, 1, MPI_LONG, MPI_SUM, communicator) ;
303
304
305  int commSize ;
306  MPI_Comm_size(communicator,&commSize) ;
307 
308        // make multiple of two
309        nbSampleNodes /= 2;
310        nbSampleNodes *= 2;
311//  assert( nbTot > nbSampleNodes*commSize) ;
312   
313  int nbSampleNodes1 = nbSampleNodes * (nb1*commSize)/(1.*nbTot) ;
314  int nbSampleNodes2 = nbSampleNodes * (nb2*commSize)/(1.*nbTot) ;
315 
316
317//      assert(node.size() > nbSampleNodes);
318//      assert(node2.size() > nbSampleNodes);
319//      assert(node.size() + node2.size() > nbSampleNodes);
320        vector<Node> sampleNodes; sampleNodes.reserve(nbSampleNodes1+nbSampleNodes2);
321
322        vector<int> randomArray1(node.size());
323        randomizeArray(randomArray1);
324        vector<int> randomArray2(node2.size());
325        randomizeArray(randomArray2);
326
327        for (int i = 0; i <nbSampleNodes1; i++) sampleNodes.push_back(Node(node[randomArray1[i%nb1]].centre,  node[randomArray1[i%nb1]].radius, NULL));
328        for (int i = 0; i <nbSampleNodes2; i++) sampleNodes.push_back(Node(node2[randomArray2[i%nb2]].centre, node2[randomArray2[i%nb2]].radius, NULL));
329
330        CTimer::get("buildParallelSampleTree").resume();
331        //sampleTree.buildParallelSampleTree(sampleNodes, cascade);
332        buildSampleTreeCascade(sampleNodes);
333        CTimer::get("buildParallelSampleTree").suspend();
334        CTimer::get("buildParallelSampleTree").print();
335
336        //route source mesh
337        CTimer::get("parallelRouteNode").resume();
338        vector<int> route(node.size());
339        cout<<"node.size = "<<node.size()<<endl;
340        routeNodes(route /*out*/, node);
341        CTimer::get("parallelRouteNode").suspend();
342        CTimer::get("parallelRouteNode").print();
343
344        CTimer::get("buildLocalTree").resume();
345        buildLocalTree(node, route);
346        CTimer::get("buildLocalTree").suspend();
347        CTimer::get("buildLocalTree").print();
348
349        CTimer::get("buildRouteTree").resume();
350        /* update circles of tree cascade so it can be used for routing */
351        updateCirclesForRouting(localTree.root->centre, localTree.root->radius);
352        CTimer::get("buildRouteTree").suspend();
353        CTimer::get("buildRouteTree").print();
354}
355
356void CParallelTree::routeNodes(vector<int>& route, vector<Node>& nodes /*route field used*/, int level)
357{
358        treeCascade[level].routeNodes(route /*assign*/, nodes, assignLevel);
359
360        if (level+1 < cascade.num_levels)
361        {
362                vector<Node> routedNodes;
363                CMPIRouting MPIRoute(cascade.level[level].pg_comm);
364                transferRoutedNodes(MPIRoute, nodes, route /*use*/, routedNodes);
365                vector<int> globalRank(routedNodes.size());
366                routeNodes(globalRank, routedNodes, level + 1);
367                MPIRoute.transferFromSource(&route[0] /*override*/, &globalRank[0]);
368        }
369        else
370        {
371                CMPIRouting MPIRoute(cascade.level[level].comm); // or use pg_comm, on last cascade level identical
372                MPIRoute.init(route);
373                int nbRecvNode = MPIRoute.getTotalSourceElement();
374                vector<int> globalRank(nbRecvNode);
375                for (int i = 0; i < globalRank.size(); i++)
376                        globalRank[i] = cascade.level[0].rank;
377                MPIRoute.transferFromSource(&route[0] /*override*/, &globalRank[0]);
378        }
379}
380
381/* assume `to` to be empty vector at entry */
382void linearize(const vector<vector<int> >& from, vector<int>& to)
383{
384        int cnt = 0;
385        for (int i = 0; i < from.size(); ++i)
386                cnt += from[i].size();
387        to.resize(cnt);
388        vector<int>::iterator dest = to.begin();
389        for (int i = 0; i < from.size(); ++i)
390                dest = copy(from[i].begin(), from[i].end(), dest);
391}
392
393/* at entry `to` must already have it's final shape and only values are overriden */
394void delinearize(const vector<int>& from, vector<vector<int> >& to)
395{
396        vector<int>::const_iterator end, src = from.begin();
397        for (int i = 0; i < to.size(); ++i, src=end)
398                copy(src, end = src + to[i].size(), to[i].begin());
399}
400
401void CParallelTree::routeIntersections(vector<vector<int> >& routes, vector<Node>& nodes, int level)
402{
403        treeCascade[level].routeIntersections(routes /*assign*/, nodes);
404
405        if (level+1 < cascade.num_levels)
406        {
407                vector<Node> routedNodes;
408                CMPIRouting MPIRoute(cascade.level[level].pg_comm);
409
410                vector<int> flattenedRoutes1;
411                linearize(routes, flattenedRoutes1);
412                vector<Node> double_nodes(flattenedRoutes1.size());
413                int j = 0;
414                for (int i = 0; i < routes.size(); ++i)
415                        for (int k = 0; k < routes[i].size(); ++k, ++j)
416                                double_nodes[j] = nodes[i];
417                transferRoutedNodes(MPIRoute, double_nodes, flattenedRoutes1 /*use*/, routedNodes);
418                vector<vector<int> > globalRanks(routedNodes.size());
419                routeIntersections(globalRanks /*out*/, routedNodes /*in*/, level + 1);
420                vector<vector<int> > flattenedRoutes(flattenedRoutes1.size());
421                // transferFromSource expects sizes (nbSendNode=flattenedRoutes, nbRecvNode=routedNodes.size())
422                MPIRoute.transferFromSource(&flattenedRoutes[0], &globalRanks[0], packVector, unpackVector);
423                for (int i = 0, j = 0; i < routes.size(); ++i)
424                {
425                        int old_size = routes[i].size();
426                        routes[i].resize(0);
427                        for (int k = 0; k < old_size; ++k, ++j)
428                                for (int l = 0; l < flattenedRoutes[j].size(); ++l)
429                                        routes[i].push_back(flattenedRoutes[j][l]);
430                }
431                assert(j == flattenedRoutes1.size());
432
433        }
434        else
435        {
436                CMPIRouting MPIRoute(cascade.level[level].comm);
437                MPIRoute.init(routes);
438                int nbRecvNode = MPIRoute.getTotalSourceElement();
439                vector<int> globalRanks(nbRecvNode, cascade.level[0].rank);
440                vector<int> flattened_routes;
441                linearize(routes, flattened_routes);
442                MPIRoute.transferFromSource(&flattened_routes[0], &globalRanks[0]);
443                delinearize(flattened_routes, routes);
444        }
445        if (level!=level)
446        {
447                for (int i = 0; i < routes.size(); ++i)
448                        for (int k = 0; k < routes[i].size(); ++k)
449                                if (routes[i][k] == cascade.level[0].rank) routes[i].erase(routes[i].begin() + (k--));
450        }
451}
452
453void CParallelTree::updateCirclesForRouting(Coord rootCentre, double rootRadius, int level)
454{
455        if (level + 1 < cascade.num_levels) // children in the cascade have to update first
456        {
457                updateCirclesForRouting(rootCentre, rootRadius, level + 1);
458                rootCentre = treeCascade[level+1].root->centre;
459                rootRadius = treeCascade[level+1].root->radius;
460        }
461
462        // gather circles on this level of the cascade
463        int pg_size;
464        MPI_Comm_size(cascade.level[level].pg_comm, &pg_size);
465        vector<Coord> allRootCentres(pg_size);
466        vector<double> allRootRadia(pg_size);
467        MPI_Allgather(&rootCentre, 3, MPI_DOUBLE, &allRootCentres[0], 3, MPI_DOUBLE, cascade.level[level].pg_comm);
468        MPI_Allgather(&rootRadius, 1, MPI_DOUBLE, &allRootRadia[0],   1, MPI_DOUBLE, cascade.level[level].pg_comm);
469
470        // now allRootsRadia and allRootCentres must be inserted into second levels of us and propagated to root
471        treeCascade[level].root->assignCircleAndPropagateUp(&allRootCentres[0], &allRootRadia[0], assignLevel);
472}
473
474CParallelTree::~CParallelTree()
475{
476        delete [] localElements;
477}
478
479}
Note: See TracBrowser for help on using the repository browser.