source: XIOS/trunk/src/transformation/transformation_mapping.cpp @ 831

Last change on this file since 831 was 831, checked in by mhnguyen, 8 years ago

Cleaning up some redundant coeds and making some improvements

+) Remove some XIOS Search to make code run faster
+) Remove some commented codes

Test
+) On Curie
+) All tests pass

File size: 9.8 KB
Line 
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#include "client_client_dht_template.hpp"
16#include "dht_data_types.hpp"
17#include "mpi_tag.hpp"
18
19namespace xios {
20
21CTransformationMapping::CTransformationMapping(CGrid* destination, CGrid* source)
22  : gridSource_(source), gridDestination_(destination)
23{
24  CContext* context = CContext::getCurrent();
25  CContextClient* client=context->client;
26  int clientRank = client->clientRank;
27
28  CDistributionClient distributionClientSrc(client->clientRank, gridSource_);
29  const CDistributionClient::GlobalLocalDataMap& globalLocalIndexGridSrc = distributionClientSrc.getGlobalLocalDataSendToServer();
30  CDistributionClient::GlobalLocalDataMap::const_iterator itIndex = globalLocalIndexGridSrc.begin(), iteIndex = globalLocalIndexGridSrc.end();
31
32  // Mapping of global index and pair containing rank and local index
33  CClientClientDHTPairIntInt::Index2InfoTypeMap globalIndexOfServer;
34  PairIntInt pairIntInt;
35  for (; itIndex != iteIndex; ++itIndex)
36  {
37    pairIntInt.first  = clientRank;
38    pairIntInt.second = itIndex->second;
39    globalIndexOfServer[itIndex->first] = pairIntInt;
40  }
41
42  gridIndexClientClientMapping_ = new CClientClientDHTPairIntInt(globalIndexOfServer,
43                                                                 client->intraComm,
44                                                                 true);
45}
46
47CTransformationMapping::CTransformationMapping(CAxis* destination, CAxis* source)
48  : gridSource_(0), gridDestination_(0)
49{
50  CContext* context = CContext::getCurrent();
51  CContextClient* client=context->client;
52  int clientRank = client->clientRank;
53
54  int niSrc     = source->n.getValue();
55  int ibeginSrc = source->begin.getValue();
56
57  CClientClientDHTPairIntInt::Index2InfoTypeMap globalIndexOfAxisSource;
58  PairIntInt pii;
59  for (int idx = 0; idx < niSrc; ++idx)
60  {
61    pii.first  = clientRank;
62    pii.second = idx;
63    globalIndexOfAxisSource[idx+ibeginSrc] = pii; //std::make_pair(clientRank,idx);
64  }
65
66  gridIndexClientClientMapping_ = new CClientClientDHTPairIntInt(globalIndexOfAxisSource,
67                                                                 client->intraComm,
68                                                                 true);
69}
70
71CTransformationMapping::~CTransformationMapping()
72{
73  if (0 != gridIndexClientClientMapping_) delete gridIndexClientClientMapping_;
74}
75
76/*!
77  Suppose that we have transformations between two grids, which are represented in form of mapping between global indexes of these two grids,
78this function tries to find out which clients a client needs to send and receive these global indexes to accomplish the transformations.
79  The grid destination is the grid whose global indexes demande global indexes from the other grid
80  Grid destination and grid source are also distributed among clients but in different manners.
81  \param [in] globaIndexWeightFromDestToSource mapping representing the transformations
82*/
83void CTransformationMapping::computeTransformationMapping(const DestinationIndexMap& globaIndexWeightFromDestToSource)
84{
85  CContext* context = CContext::getCurrent();
86  CContextClient* client=context->client;
87
88  DestinationIndexMap::const_iterator itbMap = globaIndexWeightFromDestToSource.begin(), itMap,
89                                      iteMap = globaIndexWeightFromDestToSource.end();
90
91  // Not only one index on grid destination can demande two indexes from grid source
92  // but an index on grid source has to be sent to two indexes of grid destination
93  DestinationIndexMap globalIndexMapFromSrcToDest;
94  boost::unordered_map<size_t, int> nbGlobalIndexMapFromSrcToDest;
95  std::vector<std::pair<int, std::pair<size_t,double> > >::const_iterator itbVecPair, itVecPair, iteVecPair;
96  for (itMap = itbMap; itMap != iteMap; ++itMap)
97  {
98    itbVecPair = (itMap->second).begin();
99    iteVecPair = (itMap->second).end();
100    for (itVecPair = itbVecPair; itVecPair != iteVecPair; ++itVecPair)
101    {
102      ++nbGlobalIndexMapFromSrcToDest[(itVecPair->second).first];
103    }
104  }
105
106  for (boost::unordered_map<size_t, int>::const_iterator it = nbGlobalIndexMapFromSrcToDest.begin();
107                                                         it != nbGlobalIndexMapFromSrcToDest.end(); ++it)
108  {
109    globalIndexMapFromSrcToDest[it->first].reserve(it->second);
110  }
111
112  for (itMap = itbMap; itMap != iteMap; ++itMap)
113  {
114    itbVecPair = (itMap->second).begin();
115    iteVecPair = (itMap->second).end();
116    for (itVecPair = itbVecPair; itVecPair != iteVecPair; ++itVecPair)
117    {
118      globalIndexMapFromSrcToDest[(itVecPair->second).first].push_back(std::make_pair(itVecPair->first, std::make_pair(itMap->first, (itVecPair->second).second)));
119    }
120  }
121
122  // All global indexes of a client on grid destination
123  CArray<size_t,1> globalIndexMap(globalIndexMapFromSrcToDest.size());
124  DestinationIndexMap::const_iterator itbMapBoost = globalIndexMapFromSrcToDest.begin(), itMapBoost;
125  DestinationIndexMap::const_iterator iteMapBoost = globalIndexMapFromSrcToDest.end();
126  int idx = 0;
127  for (itMapBoost = itbMapBoost; itMapBoost != iteMapBoost; ++itMapBoost)
128  {
129    globalIndexMap(idx) = itMapBoost->first;
130    ++idx;
131  }
132
133  // Find out on which clients the necessary indexes of grid source are.
134  gridIndexClientClientMapping_->computeIndexInfoMapping(globalIndexMap);
135  const CClientClientDHTPairIntInt::Index2InfoTypeMap& globalIndexSentFromGridSource = gridIndexClientClientMapping_->getInfoIndexMap();
136  CClientClientDHTPairIntInt::Index2InfoTypeMap::const_iterator itbMapSrc = globalIndexSentFromGridSource.begin(), itMapSrc,
137                                                                iteMapSrc = globalIndexSentFromGridSource.end();
138  std::vector<size_t>::const_iterator itbVec, itVec, iteVec;
139    // Inform client about the destination to which it needs to send global indexes
140  int nbClient = client->clientSize;
141  std::vector<int> sendNbClientBuff(nbClient,0);
142  std::vector<int> recvNbClientBuff(nbClient,0);
143  std::vector<int> sendIndexBuff(nbClient,0);
144  std::vector<int> recvIndexBuff(nbClient,0);
145  boost::unordered_map<int,std::vector<size_t> > sendIndexMap;
146  for (itMapSrc = itbMapSrc; itMapSrc != iteMapSrc; ++itMapSrc)
147  {
148    int sourceRank = (itMapSrc->second).first;
149    (globalIndexReceivedOnGridDestMapping_[sourceRank]).push_back(globalIndexMapFromSrcToDest[itMapSrc->first]);
150    sendIndexMap[sourceRank].push_back((itMapSrc->second).second);
151    sendIndexMap[sourceRank].push_back(itMapSrc->first);
152    sendNbClientBuff[sourceRank] = 1;
153    ++sendIndexBuff[sourceRank];
154  }
155
156  MPI_Allreduce(&sendNbClientBuff[0], &recvNbClientBuff[0], nbClient, MPI_INT, MPI_SUM, client->intraComm);
157  int numClientToReceive = recvNbClientBuff[client->clientRank];
158
159  // Then specify the size of receiving buffer, because we use synch send/receive so only necessary to know maximum size
160  MPI_Allreduce(&sendIndexBuff[0], &recvIndexBuff[0], nbClient, MPI_INT, MPI_MAX, client->intraComm);
161  int buffSize = 2*recvIndexBuff[client->clientRank]; // we send global as well as local index
162  unsigned long* recvBuffGlobalIndex;
163  if (0 != buffSize) recvBuffGlobalIndex = new unsigned long [buffSize];
164
165  std::map<int, MPI_Request> requests;
166
167  // Inform all "source clients" about index that they need to send
168  boost::unordered_map<int,std::vector<size_t> >::const_iterator itSendIndex = sendIndexMap.begin(),
169                                                                 iteSendIndex= sendIndexMap.end();
170  for (; itSendIndex != iteSendIndex; ++itSendIndex)
171  {
172    unsigned long* sendPtr = const_cast<unsigned long*>(&(itSendIndex->second)[0]);
173    MPI_Isend(sendPtr,
174              (itSendIndex->second).size(),
175              MPI_UNSIGNED_LONG,
176              (itSendIndex->first),
177              MPI_TRANSFORMATION_MAPPING_INDEX,
178              client->intraComm,
179              &requests[(itSendIndex->first)]);
180  }
181
182  // Now all the "source clients" try listening messages from other "destination clients"
183  int numClientReceived = 0;  // number of client to which data has been already sent
184  int countBuff;
185  while (numClientReceived < numClientToReceive)
186  {
187    MPI_Status status;
188    MPI_Recv(recvBuffGlobalIndex,
189             buffSize,
190             MPI_UNSIGNED_LONG,
191             MPI_ANY_SOURCE,
192             MPI_TRANSFORMATION_MAPPING_INDEX,
193             client->intraComm,
194             &status);
195
196    MPI_Get_count(&status, MPI_UNSIGNED_LONG, &countBuff);
197    int clientDestRank = status.MPI_SOURCE;
198    for (int idx = 0; idx < countBuff; idx += 2)
199    {
200      globalIndexSendToGridDestMapping_[clientDestRank].push_back(std::make_pair<int,size_t>(recvBuffGlobalIndex[idx], recvBuffGlobalIndex[idx+1]));
201    }
202    ++numClientReceived;
203  }
204
205  std::map<int, MPI_Request>::iterator itRequest;
206  for (itRequest = requests.begin(); itRequest != requests.end(); ++itRequest)
207    MPI_Wait(&itRequest->second, MPI_STATUS_IGNORE);
208
209  if (0 != buffSize) delete [] recvBuffGlobalIndex;
210}
211
212/*!
213  Return (grid) global index on grid destination. This mapping contains the rank of client source (that sends info to grid destination)
214and the corresponding global index to write on grid destination.
215  \return global index mapping to receive on grid destination
216*/
217const CTransformationMapping::ReceivedIndexMap& CTransformationMapping::getGlobalIndexReceivedOnGridDestMapping() const
218{
219  return globalIndexReceivedOnGridDestMapping_;
220}
221
222/*!
223  Return (grid) global index on grid source. This mapping contains the rank of client destination (which receives transformation info) and
224the corresponding global index to send
225  \return global index mapping to send on grid source
226*/
227const CTransformationMapping::SentIndexMap& CTransformationMapping::getGlobalIndexSendToGridDestMapping() const
228{
229  return globalIndexSendToGridDestMapping_;
230}
231
232}
Note: See TracBrowser for help on using the repository browser.