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