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

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

unify MPI_Comm type

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