source: XIOS/dev/dev_ym/XIOS_COUPLING/src/distribution/grid_remote_connector.cpp @ 2247

Last change on this file since 2247 was 2236, checked in by ymipsl, 3 years ago

Fix problem in remoteConnector when computing grid to sent to server.
Some optimisations when grid is not distributed need knowledge of the workflow view.
New CGridClientServerConnector class created based on CGridRemoteConnector.

YM

  • Property svn:eol-style set to native
  • Property svn:executable set to *
File size: 25.7 KB
Line 
1#include "grid_remote_connector.hpp"
2#include "client_client_dht_template.hpp"
3#include "leader_process.hpp"
4#include "mpi.hpp"
5
6
7
8namespace xios
9{
10  /**
11   * \brief class constructor.
12   * \param srcView List of sources views.
13   * \param dstView List of remotes views.
14   * \param localComm Local MPI communicator
15   * \param remoteSize Size of the remote communicator
16   */ 
17  CGridRemoteConnector::CGridRemoteConnector(vector<CLocalView*>& srcView, vector<CDistributedView*>& dstView, MPI_Comm localComm, int remoteSize) 
18                       : srcView_(srcView), dstView_(dstView), localComm_(localComm), remoteSize_(remoteSize) 
19  {}
20
21  /**
22   * \brief class constructor.
23   * \param srcView List of sources views.
24   * \param dstView List of remotes views.
25   * \param localComm Local MPI communicator
26   * \param remoteSize Size of the remote communicator
27   */ 
28  CGridRemoteConnector::CGridRemoteConnector(vector<CLocalView*>& srcView, vector<CLocalView*>& dstView, MPI_Comm localComm, int remoteSize) 
29                       : srcView_(srcView), localComm_(localComm), remoteSize_(remoteSize) 
30  {
31    for(auto& it : dstView) dstView_.push_back((CDistributedView*) it) ; 
32  }
33
34
35  /**
36   * \brief Compute if each view composing the source grid and the remote grid is distributed or not.
37   *         Result is stored on internal attributes \b isSrcViewDistributed_ and \b isDstViewDistributed_.
38   * \detail To compute this, a hash is computed for each array on indices. The hash must permutable, i.e.
39   *         the order of the list of global indices doesn't influence the value of the hash. So simply a sum of
40   *         hash of each indices is used for the whole array. After, the computed hash are compared with each other
41   *         ranks of \b localComm_ MPI communicator using an MPI_ALLReduce. If, for each ranks, the hash is the same
42   *         then the view is not distributed
43   */
44  void CGridRemoteConnector::computeViewDistribution(void)
45  {
46    HashXIOS<size_t> hashGlobalIndex; // hash function-object
47
48    int nDst = dstView_.size() ;
49    vector<size_t> hashRank(remoteSize_) ;
50    isDstViewDistributed_.resize(nDst) ;
51
52    for(int i=0; i<nDst; i++)
53    {
54      map<int,CArray<size_t,1>> globalIndexView ;
55      dstView_[i]->getGlobalIndexView(globalIndexView) ;
56      hashRank.assign(remoteSize_,0) ; // everybody ranks to 0 except rank of the remote view I have
57                                       // that would be assign to my local hash
58      for(auto& it : globalIndexView)
59      {
60        int rank=it.first ;
61        CArray<size_t,1>& globalIndex = it.second ;
62        size_t globalIndexSize = globalIndex.numElements();
63        size_t hashValue=0 ;
64        for(size_t ind=0;ind<globalIndexSize;ind++) hashValue += hashGlobalIndex(globalIndex(ind)) ;
65        hashRank[rank] += hashValue ;
66      }
67      // sum all the hash for every process of the local comm. The reduce is on the size of remote view (remoteSize_)
68      // after that for each rank of the remote view, we get the hash
69      MPI_Allreduce(MPI_IN_PLACE, hashRank.data(), remoteSize_, MPI_SIZE_T, MPI_SUM, localComm_) ;
70      size_t value = hashRank[0] ;
71      isDstViewDistributed_[i]=false ;
72      for(int j=0 ; j<remoteSize_ ; j++) 
73        if (value != hashRank[j]) 
74        { 
75          isDstViewDistributed_[i]=true ;
76          break ;
77        }
78    }
79
80    int nSrc = srcView_.size() ;
81    int commSize,commRank ;
82    MPI_Comm_size(localComm_,&commSize) ;
83    MPI_Comm_rank(localComm_,&commRank) ;
84    hashRank.resize(commSize,0) ;
85    isSrcViewDistributed_.resize(nSrc) ;
86
87    for(int i=0; i<nSrc; i++)
88    {
89      CArray<size_t,1> globalIndex ;
90      srcView_[i]->getGlobalIndexView(globalIndex) ;
91      hashRank.assign(commSize,0) ; // 0 for everybody except my rank
92      size_t globalIndexSize = globalIndex.numElements() ;
93      size_t hashValue=0 ;
94      for(size_t ind=0;ind<globalIndexSize;ind++) hashValue += hashGlobalIndex(globalIndex(ind)) ;
95        hashRank[commRank] += hashValue ;
96   
97      // Same method than for remote view
98      MPI_Allreduce(MPI_IN_PLACE, hashRank.data(), commSize, MPI_SIZE_T, MPI_SUM, localComm_) ;
99      size_t value = hashRank[0] ;
100      isSrcViewDistributed_[i]=false ;
101      for(int j=0 ; j<commSize ; j++) 
102        if (value != hashRank[j]) 
103        { 
104          isSrcViewDistributed_[i]=true ;
105          break ;
106        }
107    }
108
109  }
110
111/**
112  * \brief Compute the connector, i.e. compute the \b elements_ attribute.
113  * \detail Depending of the distributions of the view computed in the computeViewDistribution() call, the connector is computed in computeConnectorMethods(), and to achieve better optimisation
114  *         some redondant ranks can be removed from the elements_ map.
115  */
116  void CGridRemoteConnector::computeConnector(void)
117  {
118    computeViewDistribution() ;
119    computeConnectorMethods() ;
120    computeRedondantRanks() ; 
121    for(auto& rank : rankToRemove_)
122      for(auto& element : elements_) element.erase(rank) ;
123  }
124/**
125  * \brief Compute the connector, i.e. compute the \b elements_ attribute.
126  * \detail In order to achive better optimisation,
127  *         we distingute the case when the grid is not distributed on source grid (\bcomputeSrcNonDistributed),
128  *         or the remote grid (\b computeDstNonDistributed), or the both (\b computeSrcDstNonDistributed).
129  *         Otherwise the generic method is called computeGenericMethod. Note that in the case, if one element view
130  *         is not distributed on the source and on the remote grid, then we can used the tensorial product
131  *         property to computing it independently using \b computeSrcDstNonDistributed method.
132  *         After that, we call the \b removeRedondantRanks method to supress blocks of data that can be sent
133  *         redondantly the the remote servers
134  */
135  void CGridRemoteConnector::computeConnectorMethods(void)
136  {
137    vector<CLocalView*> srcView ;
138    vector<CDistributedView*> dstView ;
139    vector<int> indElements ;
140    elements_.resize(srcView_.size()) ;
141   
142    bool srcViewsNonDistributed=true ;
143    for(int i=0;i<srcView_.size();i++) srcViewsNonDistributed = srcViewsNonDistributed && !isSrcViewDistributed_[i]  ;
144   
145    bool dstViewsNonDistributed=true ;
146    for(int i=0;i<dstView_.size();i++) dstViewsNonDistributed = dstViewsNonDistributed && !isDstViewDistributed_[i] ;
147   
148    if (srcViewsNonDistributed) 
149    {
150      int commRank, commSize ;
151      MPI_Comm_rank(localComm_,&commRank) ;
152      MPI_Comm_size(localComm_,&commSize) ;
153      list<int> remoteRanks;
154      list<int> notUsed ;
155      map<int,bool> ranks ; 
156      computeLeaderProcess(commRank, commSize, remoteSize_, remoteRanks, notUsed) ;
157      for(int rank : remoteRanks) ranks[rank]=true ;
158     
159      for(int i=0; i<srcView_.size(); i++) 
160      {
161        if (isDstViewDistributed_[i]) computeSrcNonDistributed(i) ;
162        else computeSrcDstNonDistributed(i, ranks) ;
163      }
164    } 
165    else if (dstViewsNonDistributed)
166    {
167      map<int,bool> ranks ;
168      for(int i=0;i<remoteSize_;i++) ranks[i]=true ;
169      for(int i=0; i<srcView_.size(); i++) 
170      {
171        if (isSrcViewDistributed_[i]) computeDstNonDistributed(i,ranks) ;
172        else computeSrcDstNonDistributed(i,ranks) ;
173      }
174    } 
175    else
176    {
177      for(int i=0;i<srcView_.size();i++) 
178        if (isSrcViewDistributed_[i] || isDstViewDistributed_[i])
179        {
180          srcView.push_back(srcView_[i]) ;
181          dstView.push_back(dstView_[i]) ;
182          indElements.push_back(i) ;
183        }
184
185      computeGenericMethod(srcView, dstView, indElements) ;
186   
187      map<int,bool> ranks ; 
188      for(auto& it : elements_[indElements[0]]) 
189      {
190        if (it.second.numElements()==0) ranks[it.first] = false ;
191        else  ranks[it.first] = true ;
192      }
193   
194      for(int i=0;i<srcView_.size();i++) 
195        if (!isSrcViewDistributed_[i] && !isDstViewDistributed_[i]) computeSrcDstNonDistributed(i, ranks) ;
196    }
197
198  }
199
200 
201/**
202  * \brief Compute the connector for the element \b i when the source view is not distributed.
203  *        After the call element_[i] is defined.
204  *  \param i Indice of the element composing the source grid.
205  */
206
207  void CGridRemoteConnector::computeSrcNonDistributed(int i)
208  {
209    auto& element = elements_[i] ;
210    map<int,CArray<size_t,1>> globalIndexView ;
211    dstView_[i]->getGlobalIndexView(globalIndexView) ;
212   
213    CClientClientDHTTemplate<int>::Index2InfoTypeMap dataInfo;
214   
215    for(auto& it : globalIndexView)
216    {
217      auto& globalIndex=it.second ;
218      for(size_t ind : globalIndex) dataInfo[ind]=it.first ;
219    }
220   
221    // First we feed the distributed hash map  with key (remote global index)
222    // associated with the value of the remote rank
223    CClientClientDHTTemplate<int> DHT(dataInfo, localComm_) ;
224    // after we feed the DHT with the local global indices of the source view
225
226    int commRank, commSize ;
227    MPI_Comm_rank(localComm_,&commRank) ;
228    MPI_Comm_size(localComm_,&commSize) ;
229    CArray<size_t,1> srcIndex ;
230    // like the source view is not distributed, then only the rank 0 need to feed the DHT
231    if (commRank==0) srcView_[i]->getGlobalIndexView(srcIndex) ;
232   
233    // compute the mapping
234    DHT.computeIndexInfoMapping(srcIndex) ;
235    auto& returnInfo = DHT.getInfoIndexMap() ;
236   
237    // returnInfo contains now the map for each global indices to send to a list of remote rank
238    // only for the rank=0 because it is the one to feed the DHT
239    // so it need to send the list to each server leader i.e. the local process that handle specifically one or more
240    // servers
241   
242    // rankIndGlo : rankIndGlo[rank][indGlo] : list of indice to send the the remote server of rank "rank"
243    vector<vector<size_t>> rankIndGlo(remoteSize_) ;
244    if (commRank==0) 
245      for(auto& it1 : returnInfo)
246        for(auto& it2 : it1.second) rankIndGlo[it2].push_back(it1.first) ;
247   
248   
249    vector<MPI_Request> requests ;
250   
251    if (commRank==0)
252    {
253      requests.resize(remoteSize_) ;
254      for(int i=0 ; i<remoteSize_;i++) 
255      {
256        // ok send only the global indices for a server to the "server leader"
257        int rank = getLeaderRank(commSize, remoteSize_, i) ;
258        MPI_Isend(rankIndGlo[i].data(), rankIndGlo[i].size(), MPI_SIZE_T, rank, i ,localComm_, &requests[i]) ;
259      }
260    } 
261   
262    list<int> remoteRanks;
263    list<int> notUsed ;
264    // I am a server leader of which remote ranks ?
265    computeLeaderProcess(commRank, commSize, remoteSize_, remoteRanks, notUsed) ;
266
267    for(auto remoteRank : remoteRanks)
268    {
269      MPI_Status status ;
270      int size ;
271      MPI_Probe(0,remoteRank,localComm_, &status);
272      MPI_Get_count(&status, MPI_SIZE_T, &size) ;
273      elements_[i][remoteRank].resize(size) ;
274      // for each remote ranks receive the global indices from proc 0
275      MPI_Recv(elements_[i][remoteRank].dataFirst(),size, MPI_SIZE_T,0,remoteRank, localComm_,&status) ;
276    }
277     
278    if (commRank==0)
279    {
280      vector<MPI_Status> status(remoteSize_) ;
281      // asynchronous for sender, wait for completion
282      MPI_Waitall(remoteSize_, requests.data(), status.data()) ;
283    }
284  }
285
286  /**
287   * \brief Compute the remote connector for the element \b i when the remote view is not distributed.
288   *        After the call,  element_[i] is defined.
289   * \param i Indice of the element composing the remote grid.
290   * \param ranks The list of rank for which the local proc is in charge to compute the connector
291   *              (if leader server for exemple). if ranks[rank] == false the corresponding elements_
292   *              is set to void array (no data to sent) just in order to notify corresponding remote server
293   *              that the call is collective with each other one 
294   */
295  void CGridRemoteConnector::computeDstNonDistributed(int i, map<int,bool>& ranks)
296  {
297    auto& element = elements_[i] ;
298    map<int,CArray<size_t,1>> globalIndexView ;
299    dstView_[i]->getGlobalIndexView(globalIndexView) ;
300   
301   
302    CClientClientDHTTemplate<int>::Index2InfoTypeMap dataInfo;
303 
304    // First we feed the distributed hash map  with key (remote global index)
305    // associated with the value of the remote rank
306    for(auto& it : globalIndexView)
307      if (it.first==0) // since the remote view is not distributed, insert only the remote rank 0
308      {
309        auto& globalIndex=it.second ;
310        for(size_t ind : globalIndex) dataInfo[ind]=0 ; // associated the the rank 0
311      }
312   
313    CClientClientDHTTemplate<int> DHT(dataInfo, localComm_) ;
314    // after we feed the DHT with the local global indices of the source view
315
316    CArray<size_t,1> srcIndex ;
317    srcView_[i]->getGlobalIndexView(srcIndex) ;
318    DHT.computeIndexInfoMapping(srcIndex) ;
319    auto& returnInfo = DHT.getInfoIndexMap() ;
320   
321    // returnInfo contains now the map for each global indices to send to a list of remote rank
322    // now construct the element_ list of global indices for each rank in my list except if the erray must be empty
323    for (auto& rank : ranks)
324    {
325      if (rank.second) // non empty array => for rank that have not any data to be received
326      {
327        int size=0 ;
328        for(auto& it : returnInfo) if (!it.second.empty()) size++ ;
329        auto& array = element[rank.first] ;
330       array.resize(size) ;
331       size=0 ;
332       for(auto& it : returnInfo) 
333         if (!it.second.empty()) 
334         {
335           array(size)=it.first ;
336           size++ ;
337         }
338      }
339      else element[rank.first] = CArray<size_t,1>(0) ;  // empty array => for rank that have not any data to be received
340    }
341  }
342
343 /**
344  * \brief Compute the remote connector for the element \b i when the source and the remote view are not distributed.
345  *        After the call, element_[i] is defined.
346  * \param i Indice of the element composing the remote grid.
347  * \param ranks The list of rank for which the local proc is in charge to compute the connector
348  *              (if leader server for exemple). if ranks[rank] == false the corresponding elements_
349  *              is set to void array (no data to sent) just in order to notify corresponding remote server
350  *              that the call is collective with each other one 
351  */
352
353  void CGridRemoteConnector::computeSrcDstNonDistributed(int i, map<int,bool>& ranks)
354  {
355    auto& element = elements_[i] ;
356    map<int,CArray<size_t,1>> globalIndexView ;
357    dstView_[i]->getGlobalIndexView(globalIndexView) ;
358   
359   
360    CClientClientDHTTemplate<int>::Index2InfoTypeMap dataInfo;
361    // First we feed the distributed hash map  with key (remote global index)
362    // associated with the value of the remote rank
363
364    for(auto& it : globalIndexView)
365      if (it.first==0) // insert only the remote rank 0 since the remote view is not distributed
366      {
367        auto& globalIndex=it.second ;
368        for(size_t ind : globalIndex) dataInfo[ind]=0 ; // associated the the rank 0
369      }
370   
371    CClientClientDHTTemplate<int> DHT(dataInfo, localComm_) ;
372    // after we feed the DHT with the local global indices of the source view
373
374    int commRank, commSize ;
375    MPI_Comm_rank(localComm_,&commRank) ;
376    MPI_Comm_size(localComm_,&commSize) ;
377    CArray<size_t,1> srcIndex ;
378 
379    // like the source view is not distributed, then only the rank 0 need to feed the DHT
380    if (commRank==0) srcView_[i]->getGlobalIndexView(srcIndex) ;
381    DHT.computeIndexInfoMapping(srcIndex) ;
382    auto& returnInfo = DHT.getInfoIndexMap() ;
383   
384    vector<size_t> indGlo ;
385    if (commRank==0) 
386      for(auto& it1 : returnInfo) 
387        for(auto& it2 : it1.second) indGlo.push_back(it1.first) ;
388
389    // now local rank 0 know which indices to seed to remote rank 0, but all the server
390    // must receive the same information. So only the leader rank will sent this.
391    // So local rank 0 must broadcast the information to all leader.
392    // for this we create a new communicator composed of local process that must send data
393    // to a remote rank, data are broadcasted, and element_[i] is construction for each remote
394    // rank in charge
395    int color=0 ;
396    if (ranks.empty()) color=0 ;
397    else color=1 ;
398    if (commRank==0) color=1 ;
399    MPI_Comm newComm ;
400    MPI_Comm_split(localComm_, color, commRank, &newComm) ;
401    if (color==1)
402    {
403      // ok, I am part of the process that must send something to one or more remote server
404      // so I get the list of global indices from rank 0
405      int dataSize ;
406      if (commRank==0) dataSize=indGlo.size() ;
407      MPI_Bcast(&dataSize,1,MPI_INT, 0, newComm) ;
408      indGlo.resize(dataSize) ;
409      MPI_Bcast(indGlo.data(),dataSize,MPI_SIZE_T,0,newComm) ;
410    }
411    MPI_Comm_free(&newComm) ;
412
413    // construct element_[i] from indGlo
414    for(auto& rank : ranks)
415    {
416      if (rank.second)
417      {
418        int dataSize=indGlo.size();
419        auto& element = elements_[i][rank.first] ;
420        element.resize(dataSize) ;
421        for(int i=0;i<dataSize; i++) element(i)=indGlo[i] ;
422      }
423      else element[rank.first] = CArray<size_t,1>(0) ;
424    }   
425
426  }
427
428 
429 /**
430  * \brief Generic method the compute the grid remote connector. Only distributed elements are specifed in the source view and remote view.
431  *        Connector for non distributed elements are computed separatly to improve performance and memory consumption. After the call,
432  *        \b elements_  is defined.
433  *  \param srcView List of the source views composing the grid, without non distributed views
434  *  \param dstView List of the remote views composing the grid, without non distributed views
435  *  \param indElements Index of the view making the correspondance between all views and views distributed (that are in input)
436  */
437  void CGridRemoteConnector::computeGenericMethod(vector<CLocalView*>& srcView, vector<CDistributedView*>& dstView, vector<int>& indElements)
438  {
439    // generic method, every element can be distributed
440    int nDst = dstView.size() ;
441    vector<size_t> dstSliceSize(nDst) ;
442    dstSliceSize[0] = 1 ; 
443    for(int i=1; i<nDst; i++)  dstSliceSize[i] = dstView[i-1]->getGlobalSize()*dstSliceSize[i-1] ;
444 
445    CClientClientDHTTemplate<int>::Index2VectorInfoTypeMap dataInfo ;
446    CClientClientDHTTemplate<size_t>::Index2VectorInfoTypeMap info ; // info map
447
448    // first, we need to feed the DHT with the global index of the remote server
449    // for that :
450    // First the first element insert the in a DHT with key as the rank and value the list of global index associated
451    // Then get the previously stored index associate with the remote rank I am in charge and reinsert the global index
452    // corresponding to the position of the element in the remote view suing tensorial product
453    // finaly we get only the list of remote global index I am in charge for the whole remote grid   
454
455    for(int pos=0; pos<nDst; pos++)
456    {
457      size_t sliceSize=dstSliceSize[pos] ;
458      map<int,CArray<size_t,1>> globalIndexView ;
459      dstView[pos]->getGlobalIndexView(globalIndexView) ;
460     
461      CClientClientDHTTemplate<size_t>::Index2VectorInfoTypeMap lastInfo(info) ;
462
463      if (pos>0)
464      {
465        CArray<size_t,1> ranks(globalIndexView.size()) ;
466        auto it=globalIndexView.begin() ;
467        for(int i=0 ; i<ranks.numElements();i++,it++) ranks(i)=it->first ;
468        CClientClientDHTTemplate<size_t> dataRanks(info, localComm_) ;
469        dataRanks.computeIndexInfoMapping(ranks) ;
470        lastInfo = dataRanks.getInfoIndexMap() ;
471      }
472     
473      info.clear() ;
474      for(auto& it : globalIndexView)
475      {
476        int rank = it.first ;
477        auto& globalIndex = it.second ;
478        auto& inf = info[rank] ;
479        if (pos==0) for(int i=0;i<globalIndex.numElements();i++) inf.push_back(globalIndex(i)) ;
480        else
481        {
482          auto& lastGlobalIndex = lastInfo[rank] ;
483          for(size_t lastGlobalInd : lastGlobalIndex)
484          {
485            for(int i=0;i<globalIndex.numElements();i++) inf.push_back(globalIndex(i)*sliceSize+lastGlobalInd) ;
486          }
487        } 
488      }
489
490      if (pos==nDst-1)
491      {
492         for(auto& it : info)
493         {
494           int rank=it.first ;
495           auto& globalIndex = it.second ;
496           for(auto globalInd : globalIndex) dataInfo[globalInd].push_back(rank) ;
497         }
498      } 
499    }
500
501    // we feed the DHT with the remote global index
502    CClientClientDHTTemplate<int> dataRanks(dataInfo, localComm_) ;
503
504    // generate list of global index for src view
505    int nSrc = srcView.size() ;
506    vector<size_t> srcSliceSize(nSrc) ;
507   
508    srcSliceSize[0] = 1 ; 
509    for(int i=1; i<nSrc; i++)  srcSliceSize[i] = srcView[i-1]->getGlobalSize()*srcSliceSize[i-1] ;
510
511    vector<size_t> srcGlobalIndex ;
512    size_t sliceIndex=0 ;
513    srcView[nSrc-1]->getGlobalIndex(srcGlobalIndex, sliceIndex, srcSliceSize.data(), srcView.data(), nSrc-1) ;
514    // now we have the global index of the source grid in srcGlobalIndex
515    // we feed the DHT with the src global index (if we have)
516    if (srcGlobalIndex.size()>0)
517    {
518      CArray<size_t,1> srcGlobalIndexArray(srcGlobalIndex.data(), shape(srcGlobalIndex.size()),neverDeleteData) ;
519      dataRanks.computeIndexInfoMapping(srcGlobalIndexArray) ;
520    }
521    else
522    {
523      CArray<size_t,1> srcGlobalIndexArray ;
524      dataRanks.computeIndexInfoMapping(srcGlobalIndexArray) ;
525    }
526    const auto& returnInfo = dataRanks.getInfoIndexMap() ;
527    // returnInfo contains now the map for each global indices to send to a list of remote rank
528    // but we want to use the tensorial product property to get the same information using only global
529    // index of element view. So the idea is to reverse the information : for a global index of the grid
530    // to send to the remote server, what is the global index of each element composing the grid ?
531
532    vector<map<int, set<size_t>>> elements(nSrc) ; // internal representation of elements composing the grid
533
534    for(auto& indRanks : returnInfo)
535    {
536      size_t gridIndexGlo=indRanks.first ;
537      auto& ranks = indRanks.second ;
538      for(int i=nSrc-1; i>=0; i--)
539      {
540        auto& element = elements[i] ;
541        size_t localIndGlo = gridIndexGlo / srcSliceSize[i] ;
542        gridIndexGlo = gridIndexGlo % srcSliceSize[i] ;
543        for(int rank : ranks) element[rank].insert(localIndGlo) ;
544      }
545    }
546
547//    elements_.resize(nSrc) ;
548    for(int i=0 ; i<nSrc; i++)
549    {
550      auto& element=elements[i] ;
551      for(auto& rankInd : element)
552      {
553        int rank=rankInd.first ;
554        set<size_t>& indGlo = rankInd.second ;
555        CArray<size_t,1>& indGloArray = elements_[indElements[i]][rank] ;
556        indGloArray.resize(indGlo.size()) ;
557        int j=0 ;
558        for (auto index : indGlo) { indGloArray(j) = index ; j++; }
559      }
560    }
561   
562    // So what about when there is some server that have no data to receive
563    // they must be inform they receive an event with no data.
564    // So find remote servers with no data, and one client will take in charge
565    // that it receive global index with no data (0-size)
566    vector<int> ranks(remoteSize_,0) ;
567    for(auto& it : elements_[indElements[0]]) ranks[it.first] = 1 ;
568    MPI_Allreduce(MPI_IN_PLACE, ranks.data(), remoteSize_, MPI_INT, MPI_SUM, localComm_) ;
569    int commRank, commSize ;
570    MPI_Comm_rank(localComm_,&commRank) ;
571    MPI_Comm_size(localComm_,&commSize) ;
572    int pos=0 ;
573    for(int i=0; i<remoteSize_ ; i++)
574      if (ranks[i]==0)
575      {
576        if (pos%commSize==commRank) 
577          for(int j=0 ; j<nSrc; j++) elements_[indElements[j]][i] = CArray<size_t,1>(0) ;
578        pos++ ;
579      }
580  }
581
582 /**
583  * \brief Once the connector is computed (compute \b elements_), redondant data can be avoid to be sent to the server.
584  *        This call compute the redondant rank and store them in \b rankToRemove_ attribute.
585  *        The goal of this method is to make a hash of each block of indice that determine wich data to send to a
586  *        of a specific server rank using a hash method. So data to send to a rank is associated to a hash.
587  *        After we compare hash between local rank and remove redondant data corresponding to the same hash.
588  */
589  void CGridRemoteConnector::computeRedondantRanks(void)
590  {
591    int commRank ;
592    MPI_Comm_rank(localComm_,&commRank) ;
593
594    set<int> ranks;
595    for(auto& element : elements_) 
596      for(auto& it : element) ranks.insert(it.first) ;
597
598    for(auto& element : elements_)
599      for(auto& it : element) 
600        if (ranks.count(it.first)==0) ERROR("void CGridRemoteConnector::removeRedondantRanks(void)",<<"number of ranks in elements is not coherent between each element") ;
601   
602    HashXIOS<size_t> hashGlobalIndex;
603   
604    map<int,size_t> hashRanks ;
605    for(auto& element : elements_) 
606      for(auto& it : element)
607      {
608        auto& globalIndex=it.second ;
609        int rank=it.first ;
610        size_t hash ;
611        hash=0 ;
612        for(int i=0; i<globalIndex.numElements(); i++) hash+=hashGlobalIndex(globalIndex(i)) ;
613        if (hashRanks.count(rank)==0) hashRanks[rank]=hash ;
614        else hashRanks[rank]=hashGlobalIndex.hashCombine(hashRanks[rank],hash) ;
615      }
616    // a hash is now computed for data block I will sent to the server.
617
618    CClientClientDHTTemplate<int>::Index2InfoTypeMap info ;
619
620    map<size_t,int> hashRank ;
621    HashXIOS<int> hashGlobalIndexRank;
622    for(auto& it : hashRanks) 
623    {
624      it.second = hashGlobalIndexRank.hashCombine(it.first,it.second) ; 
625      info[it.second]=commRank ;
626      hashRank[it.second]=it.first ;
627    }
628
629    // we feed a DHT map with key : hash, value : myrank
630    CClientClientDHTTemplate<int> dataHash(info, localComm_) ;
631    CArray<size_t,1> hashList(hashRank.size()) ;
632   
633    int i=0 ;
634    for(auto& it : hashRank) { hashList(i)=it.first ; i++; }
635
636    // now who are the ranks that have the same hash : feed the DHT with my list of hash
637    dataHash.computeIndexInfoMapping(hashList) ;
638    auto& hashRankList = dataHash.getInfoIndexMap() ;
639   
640
641    for(auto& it : hashRankList)
642    {
643      size_t hash = it.first ;
644      auto& ranks = it.second ;
645     
646      bool first=true ;
647      // only the process with the lowest rank get in charge of sendinf data to remote server
648      for(int rank : ranks) if (commRank>rank) first=false ;
649      if (!first) rankToRemove_.insert(hashRank[hash]) ;
650    }
651  }
652
653}
Note: See TracBrowser for help on using the repository browser.