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

Last change on this file since 2296 was 2296, checked in by ymipsl, 2 years ago

Solve problem rising when output scalar.
Jenkins hash method map 0 to 0 which conflict with hash initialisation that is also set to 0.

YM

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