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

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

intermediate commit for new tranformation engine?
YM

  • Property svn:eol-style set to native
  • Property svn:executable set to *
File size: 5.0 KB
Line 
1#include "grid_remote_connector.hpp"
2#include "client_client_dht_template.hpp"
3#include "mpi.hpp"
4
5
6
7namespace xios
8{
9 
10  CGridRemoteConnector::CGridRemoteConnector(vector<CLocalView*>& srcView, vector<CDistributedView*>& dstView, MPI_Comm localComm, int remoteSize) 
11                       : srcView_(srcView), dstView_(dstView), localComm_(localComm), remoteSize_(remoteSize) 
12  {}
13
14  void CGridRemoteConnector::computeConnector(void)
15  {
16    computeGenericMethod() ;
17  }
18
19  void CGridRemoteConnector::computeGenericMethod(void)
20  {
21    // generic method, every element can be distributed
22    int nDst = dstView_.size() ;
23    vector<size_t> dstSliceSize(nDst) ;
24    dstSliceSize[0] = 1 ; 
25    for(int i=1; i<nDst; i++)  dstSliceSize[i] = dstView_[i-1]->getGlobalSize()*dstSliceSize[i-1] ;
26 
27    CClientClientDHTTemplate<int>::Index2VectorInfoTypeMap dataInfo ;
28
29    CClientClientDHTTemplate<size_t>::Index2VectorInfoTypeMap info ; // info map
30    for(int pos=0; pos<nDst; pos++)
31    {
32      size_t sliceSize=dstSliceSize[pos] ;
33      map<int,CArray<size_t,1>> globalIndexView ;
34      dstView_[pos]->getGlobalIndexView(globalIndexView) ;
35     
36      CClientClientDHTTemplate<size_t>::Index2VectorInfoTypeMap lastInfo(info) ;
37
38      if (pos>0)
39      {
40        CArray<size_t,1> ranks(globalIndexView.size()) ;
41        auto it=globalIndexView.begin() ;
42        for(int i=0 ; i<ranks.numElements();i++,it++) ranks(i)=it->first ;
43        CClientClientDHTTemplate<size_t> dataRanks(info, localComm_) ;
44        dataRanks.computeIndexInfoMapping(ranks) ;
45        lastInfo = dataRanks.getInfoIndexMap() ;
46      }
47     
48      info.clear() ;
49      for(auto& it : globalIndexView)
50      {
51        int rank = it.first ;
52        auto& globalIndex = it.second ;
53        auto& inf = info[rank] ;
54        if (pos==0) for(int i=0;i<globalIndex.numElements();i++) inf.push_back(globalIndex(i)) ;
55        else
56        {
57          auto& lastGlobalIndex = lastInfo[rank] ;
58          for(size_t lastGlobalInd : lastGlobalIndex)
59          {
60            for(int i=0;i<globalIndex.numElements();i++) inf.push_back(globalIndex(i)*sliceSize+lastGlobalInd) ;
61          }
62        } 
63      }
64
65      if (pos==nDst-1)
66      {
67         for(auto& it : info)
68         {
69           int rank=it.first ;
70           auto& globalIndex = it.second ;
71           for(auto globalInd : globalIndex) dataInfo[globalInd].push_back(rank) ;
72         }
73      } 
74    }
75
76    CClientClientDHTTemplate<int> dataRanks(dataInfo, localComm_) ;
77
78    // generate list of global index for src view
79    int nSrc = srcView_.size() ;
80    vector<size_t> srcSliceSize(nSrc) ;
81   
82    srcSliceSize[0] = 1 ; 
83    for(int i=1; i<nSrc; i++)  srcSliceSize[i] = srcView_[i-1]->getGlobalSize()*srcSliceSize[i-1] ;
84
85    vector<size_t> srcGlobalIndex ;
86    size_t sliceIndex=0 ;
87    srcView_[nSrc-1]->getGlobalIndex(srcGlobalIndex, sliceIndex, srcSliceSize.data(), srcView_.data(), nSrc-1) ;
88   
89    if (srcGlobalIndex.size()>0)
90    {
91      CArray<size_t,1> srcGlobalIndexArray(srcGlobalIndex.data(), shape(srcGlobalIndex.size()),neverDeleteData) ;
92      dataRanks.computeIndexInfoMapping(srcGlobalIndexArray) ;
93    }
94    else
95    {
96      CArray<size_t,1> srcGlobalIndexArray ;
97      dataRanks.computeIndexInfoMapping(srcGlobalIndexArray) ;
98    }
99    const auto& returnInfo = dataRanks.getInfoIndexMap() ;
100
101    vector<map<int, set<size_t>>> elements(nSrc) ; // internal representation of elements composing the grid
102
103    for(auto& indRanks : returnInfo)
104    {
105      size_t gridIndexGlo=indRanks.first ;
106      auto& ranks = indRanks.second ;
107      for(int i=nSrc-1; i>=0; i--)
108      {
109        auto& element = elements[i] ;
110        size_t localIndGlo = gridIndexGlo / srcSliceSize[i] ;
111        gridIndexGlo = gridIndexGlo % srcSliceSize[i] ;
112        for(int rank : ranks) element[rank].insert(localIndGlo) ;
113      }
114    }
115
116    elements_.resize(nSrc) ;
117    for(int i=0 ; i<nSrc; i++)
118    {
119      auto& element=elements[i] ;
120      for(auto& rankInd : element)
121      {
122        int rank=rankInd.first ;
123        set<size_t>& indGlo = rankInd.second ;
124        CArray<size_t,1>& indGloArray = elements_[i][rank] ;
125        indGloArray.resize(indGlo.size()) ;
126        int j=0 ;
127        for (auto index : indGlo) { indGloArray(j) = index ; j++; }
128      }
129    }
130   
131    // So what about when there is some server that have no data to receive
132    // they must be inform they receive an event with no data.
133    // So find remote servers with no data, and one client will take in charge
134    // that it receive global index with no data (0-size)
135    vector<int> ranks(remoteSize_,0) ;
136    for(auto& it : elements_[0]) ranks[it.first] = 1 ;
137    MPI_Allreduce(MPI_IN_PLACE, ranks.data(), remoteSize_, MPI_INT, MPI_SUM, localComm_) ;
138    int commRank, commSize ;
139    MPI_Comm_rank(localComm_,&commRank) ;
140    MPI_Comm_size(localComm_,&commSize) ;
141    int pos=0 ;
142    for(int i=0; i<remoteSize_ ; i++)
143      if (ranks[i]==0)
144      {
145        if (pos%commSize==commRank) for(auto& element : elements_) element[i] = CArray<size_t,1>(0) ;
146        pos++ ;
147      }
148  }
149
150
151}
Note: See TracBrowser for help on using the repository browser.