1 | /*! |
---|
2 | \file transformation_mapping.cpp |
---|
3 | \author Ha NGUYEN |
---|
4 | \since 14 May 2015 |
---|
5 | \date 18 June 2015 |
---|
6 | |
---|
7 | \brief Take charge of communication among clients to exchange transformed data. |
---|
8 | */ |
---|
9 | |
---|
10 | #include "transformation_mapping.hpp" |
---|
11 | #include <boost/unordered_map.hpp> |
---|
12 | #include "context.hpp" |
---|
13 | #include "context_client.hpp" |
---|
14 | #include "distribution_client.hpp" |
---|
15 | |
---|
16 | namespace xios { |
---|
17 | |
---|
18 | CTransformationMapping::CTransformationMapping(CGrid* destination, CGrid* source) |
---|
19 | : gridSource_(source), gridDestination_(destination) |
---|
20 | { |
---|
21 | CContext* context = CContext::getCurrent(); |
---|
22 | CContextClient* client=context->client; |
---|
23 | int clientRank = client->clientRank; |
---|
24 | |
---|
25 | CDistributionClient distributionClientSrc(client->clientRank, gridSource_); |
---|
26 | |
---|
27 | const CArray<size_t,1>& globalIndexGridSrc = distributionClientSrc.getGlobalDataIndexSendToServer(); //gridSource_->getDistributionClient()->getGlobalDataIndexSendToServer(); |
---|
28 | boost::unordered_map<size_t,int> globalIndexOfServer; |
---|
29 | int globalIndexSize = globalIndexGridSrc.numElements(); |
---|
30 | for (int idx = 0; idx < globalIndexSize; ++idx) |
---|
31 | { |
---|
32 | globalIndexOfServer[globalIndexGridSrc(idx)] = clientRank; |
---|
33 | } |
---|
34 | |
---|
35 | gridIndexClientClientMapping_ = new CClientServerMappingDistributed(globalIndexOfServer, |
---|
36 | client->intraComm, |
---|
37 | true); |
---|
38 | } |
---|
39 | |
---|
40 | CTransformationMapping::CTransformationMapping(CAxis* destination, CAxis* source) |
---|
41 | : gridSource_(0), gridDestination_(0) |
---|
42 | { |
---|
43 | CContext* context = CContext::getCurrent(); |
---|
44 | CContextClient* client=context->client; |
---|
45 | int clientRank = client->clientRank; |
---|
46 | |
---|
47 | int niSrc = source->ni.getValue(); |
---|
48 | int ibeginSrc = source->ibegin.getValue(); |
---|
49 | |
---|
50 | boost::unordered_map<size_t,int> globalIndexOfAxisSource; |
---|
51 | for (int idx = 0; idx < niSrc; ++idx) |
---|
52 | { |
---|
53 | globalIndexOfAxisSource[idx+ibeginSrc] = clientRank; |
---|
54 | } |
---|
55 | |
---|
56 | gridIndexClientClientMapping_ = new CClientServerMappingDistributed(globalIndexOfAxisSource, |
---|
57 | client->intraComm, |
---|
58 | true); |
---|
59 | } |
---|
60 | |
---|
61 | CTransformationMapping::~CTransformationMapping() |
---|
62 | { |
---|
63 | if (0 != gridIndexClientClientMapping_) delete gridIndexClientClientMapping_; |
---|
64 | } |
---|
65 | |
---|
66 | /*! |
---|
67 | Suppose that we have transformations between two grids, which are represented in form of mapping between global indexes of these two grids, |
---|
68 | this function tries to find out which clients a client needs to send and receive these global indexes to accomplish the transformations. |
---|
69 | The grid destination is the grid whose global indexes demande global indexes from the other grid |
---|
70 | Grid destination and grid source are also distributed among clients but in different manners. |
---|
71 | \param [in] globaIndexWeightFromDestToSource mapping representing the transformations |
---|
72 | */ |
---|
73 | void CTransformationMapping::computeTransformationMapping(const std::map<size_t, std::vector<std::pair<size_t,double> > >& globaIndexWeightFromDestToSource) |
---|
74 | { |
---|
75 | CContext* context = CContext::getCurrent(); |
---|
76 | CContextClient* client=context->client; |
---|
77 | |
---|
78 | std::map<size_t, std::vector<std::pair<size_t,double> > >::const_iterator itbMap = globaIndexWeightFromDestToSource.begin(), itMap, |
---|
79 | iteMap = globaIndexWeightFromDestToSource.end(); |
---|
80 | |
---|
81 | // Not only one index on grid destination can demande two indexes from grid source |
---|
82 | // but an index on grid source has to be sent to two indexes of grid destination |
---|
83 | std::map<size_t, std::vector<std::pair<size_t,double> > > globalIndexMapFromSrcToDest; |
---|
84 | std::vector<std::pair<size_t,double> >::const_iterator itbVecPair, itVecPair, iteVecPair; |
---|
85 | for (itMap = itbMap; itMap != iteMap; ++itMap) |
---|
86 | { |
---|
87 | itbVecPair = (itMap->second).begin(); |
---|
88 | iteVecPair = (itMap->second).end(); |
---|
89 | for (itVecPair = itbVecPair; itVecPair != iteVecPair; ++itVecPair) |
---|
90 | { |
---|
91 | globalIndexMapFromSrcToDest[itVecPair->first].push_back(std::make_pair(itMap->first, itVecPair->second)); |
---|
92 | } |
---|
93 | } |
---|
94 | |
---|
95 | // All global indexes of a client on grid destination |
---|
96 | CArray<size_t,1> globalIndexMap(globalIndexMapFromSrcToDest.size()); |
---|
97 | itbMap = globalIndexMapFromSrcToDest.begin(); |
---|
98 | iteMap = globalIndexMapFromSrcToDest.end(); |
---|
99 | int idx = 0; |
---|
100 | for (itMap = itbMap; itMap != iteMap; ++itMap) |
---|
101 | { |
---|
102 | globalIndexMap(idx) = itMap->first; |
---|
103 | ++idx; |
---|
104 | } |
---|
105 | |
---|
106 | // Find out on which clients the necessary indexes of grid source are. |
---|
107 | gridIndexClientClientMapping_->computeServerIndexMapping(globalIndexMap); |
---|
108 | const std::map<int, std::vector<size_t> >& globalIndexSentFromGridSource = gridIndexClientClientMapping_->getGlobalIndexOnServer(); |
---|
109 | std::map<int, std::vector<size_t> >::const_iterator itbMapSrc = globalIndexSentFromGridSource.begin(), itMapSrc, |
---|
110 | iteMapSrc = globalIndexSentFromGridSource.end(); |
---|
111 | std::vector<size_t>::const_iterator itbVec, itVec, iteVec; |
---|
112 | for (itMapSrc = itbMapSrc; itMapSrc != iteMapSrc; ++itMapSrc) |
---|
113 | { |
---|
114 | int sourceRank = itMapSrc->first; |
---|
115 | itbVec = (itMapSrc->second).begin(); |
---|
116 | iteVec = (itMapSrc->second).end(); |
---|
117 | for (itVec = itbVec; itVec != iteVec; ++itVec) |
---|
118 | { |
---|
119 | (globalIndexReceivedOnGridDestMapping_[sourceRank]).push_back(globalIndexMapFromSrcToDest[*itVec]); |
---|
120 | } |
---|
121 | } |
---|
122 | |
---|
123 | // Inform client about the destination to which it needs to send global indexes |
---|
124 | int nbClient = client->clientSize; |
---|
125 | int* sendBuff = new int[nbClient]; |
---|
126 | int* recvBuff = new int[nbClient]; |
---|
127 | for (int i = 0; i < nbClient; ++i) sendBuff[i] = 0; |
---|
128 | |
---|
129 | // First of all, inform the number of destination a client needs to send global index |
---|
130 | for (itMapSrc = itbMapSrc; itMapSrc != iteMapSrc; ++itMapSrc) sendBuff[itMapSrc->first] = 1; |
---|
131 | MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm); |
---|
132 | int numClientToReceive = recvBuff[client->clientRank]; |
---|
133 | |
---|
134 | // Then specify the size of receiving buffer, because we use synch send/receive so only necessary to know maximum size |
---|
135 | for (itMapSrc = itbMapSrc; itMapSrc != iteMapSrc; ++itMapSrc) sendBuff[itMapSrc->first] = (itMapSrc->second).size(); |
---|
136 | MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_MAX, client->intraComm); |
---|
137 | |
---|
138 | int buffSize = recvBuff[client->clientRank]; |
---|
139 | unsigned long* recvBuffGlobalIndex; |
---|
140 | if (0 != buffSize) recvBuffGlobalIndex = new unsigned long [buffSize]; |
---|
141 | |
---|
142 | std::map<int, MPI_Request> requests; |
---|
143 | |
---|
144 | // Inform all "source clients" about index that they need to send |
---|
145 | for (itMapSrc = itbMapSrc; itMapSrc != iteMapSrc; ++itMapSrc) |
---|
146 | { |
---|
147 | unsigned long* sendPtr = const_cast<unsigned long*>(&(itMapSrc->second)[0]); |
---|
148 | MPI_Isend(sendPtr, |
---|
149 | (itMapSrc->second).size(), |
---|
150 | MPI_UNSIGNED_LONG, |
---|
151 | itMapSrc->first, |
---|
152 | 11, |
---|
153 | client->intraComm, |
---|
154 | &requests[itMapSrc->first]); |
---|
155 | } |
---|
156 | |
---|
157 | // Now all the "source clients" try listening messages from other "destination clients" |
---|
158 | int numClientReceived = 0; // number of client to which data has been already sent |
---|
159 | int countBuff; |
---|
160 | while (numClientReceived < numClientToReceive) |
---|
161 | { |
---|
162 | MPI_Status status; |
---|
163 | MPI_Recv(recvBuffGlobalIndex, |
---|
164 | buffSize, |
---|
165 | MPI_UNSIGNED_LONG, |
---|
166 | MPI_ANY_SOURCE, |
---|
167 | 11, |
---|
168 | client->intraComm, |
---|
169 | &status); |
---|
170 | |
---|
171 | MPI_Get_count(&status, MPI_UNSIGNED_LONG, &countBuff); |
---|
172 | int clientDestRank = status.MPI_SOURCE; |
---|
173 | for (int idx = 0; idx < countBuff; ++idx) |
---|
174 | { |
---|
175 | globalIndexSendToGridDestMapping_[clientDestRank].push_back(recvBuffGlobalIndex[idx]); |
---|
176 | } |
---|
177 | ++numClientReceived; |
---|
178 | } |
---|
179 | |
---|
180 | std::map<int, MPI_Request>::iterator itRequest; |
---|
181 | for (itRequest = requests.begin(); itRequest != requests.end(); ++itRequest) |
---|
182 | MPI_Wait(&itRequest->second, MPI_STATUS_IGNORE); |
---|
183 | |
---|
184 | delete [] sendBuff; |
---|
185 | delete [] recvBuff; |
---|
186 | if (0 != buffSize) delete [] recvBuffGlobalIndex; |
---|
187 | } |
---|
188 | |
---|
189 | /*! |
---|
190 | Return (grid) global index on grid destination. This mapping contains the rank of client source (that sends info to grid destination) |
---|
191 | and the corresponding global index to write on grid destination. |
---|
192 | \return global index mapping to receive on grid destination |
---|
193 | */ |
---|
194 | const std::map<int,std::vector<std::vector<std::pair<size_t,double> > > >& CTransformationMapping::getGlobalIndexReceivedOnGridDestMapping() const |
---|
195 | { |
---|
196 | return globalIndexReceivedOnGridDestMapping_; |
---|
197 | } |
---|
198 | |
---|
199 | /*! |
---|
200 | Return (grid) global index on grid source. This mapping contains the rank of client destination (which receives transformation info) and |
---|
201 | the corresponding global index to send |
---|
202 | \return global index mapping to send on grid source |
---|
203 | */ |
---|
204 | const std::map<int,std::vector<size_t> >& CTransformationMapping::getGlobalIndexSendToGridDestMapping() const |
---|
205 | { |
---|
206 | return globalIndexSendToGridDestMapping_; |
---|
207 | } |
---|
208 | |
---|
209 | } |
---|