source: XIOS/trunk/src/transformation/domain_algorithm_interpolate.cpp @ 689

Last change on this file since 689 was 689, checked in by mhnguyen, 9 years ago

Modifying the interface of interpolation domain

+) Change node name from interpolate_from_file_domain to interpolate_domain and add some new atrributes
+) Add more tests into test_remap

Test
+) On Curie
+) test_remap works for direct weight calculation and reading weight calculation from file

File size: 12.1 KB
Line 
1/*!
2   \file domain_algorithm_interpolate_from_file.cpp
3   \author Ha NGUYEN
4   \since 09 Jul 2015
5   \date 15 Sep 2015
6
7   \brief Algorithm for interpolation on a domain.
8 */
9#include "domain_algorithm_interpolate.hpp"
10#include <boost/unordered_map.hpp>
11#include "context.hpp"
12#include "context_client.hpp"
13#include "distribution_client.hpp"
14#include "client_server_mapping_distributed.hpp"
15#include "netcdf.hpp"
16#include "mapper.hpp"
17
18namespace xios {
19
20CDomainAlgorithmInterpolate::CDomainAlgorithmInterpolate(CDomain* domainDestination, CDomain* domainSource, CInterpolateDomain* interpDomain)
21: CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain)
22{
23  interpDomain_->checkValid(domainSource);
24  computeIndexSourceMapping();
25}
26
27/*!
28  Compute remap with integrated remap calculation module
29*/
30void CDomainAlgorithmInterpolate::computeRemap()
31{
32  using namespace sphereRemap;
33
34  CContext* context = CContext::getCurrent();
35  CContextClient* client=context->client;
36  int clientRank = client->clientRank;
37  int i, j, k, idx;
38  std::vector<double> srcPole(3,0), dstPole(3,0);
39        int orderInterp = interpDomain_->order.getValue();
40
41  int constNVertex = 4; // Value by default number of vertex for rectangular domain
42  int nVertexSrc, nVertexDest;
43  nVertexSrc = nVertexDest = constNVertex;
44
45  // First of all, try to retrieve the boundary values of domain source and domain destination
46  int localDomainSrcSize = domainSrc_->i_index.numElements();
47  int niSrc = domainSrc_->ni.getValue(), njSrc = domainSrc_->nj.getValue();
48  bool hasBoundSrc = domainSrc_->hasBounds;
49  if (hasBoundSrc) nVertexSrc = domainSrc_->nvertex.getValue();
50  CArray<double,2> boundsLonSrc(nVertexSrc,localDomainSrcSize);
51  CArray<double,2> boundsLatSrc(nVertexSrc,localDomainSrcSize);
52
53  if (CDomain::type_attr::rectilinear == domainSrc_->type) srcPole[2] = 1;
54  if (hasBoundSrc)  // Suppose that domain source is curvilinear or unstructured
55  {
56    if (!domainSrc_->bounds_lon_2d.isEmpty())
57    {
58      for (j = 0; j < njSrc; ++j)
59        for (i = 0; i < niSrc; ++i)
60        {
61          k=j*niSrc+i;
62          for(int n=0;n<nVertexSrc;++n)
63          {
64            boundsLonSrc(n,k) = domainSrc_->bounds_lon_2d(n,i,j);
65            boundsLatSrc(n,k) = domainSrc_->bounds_lat_2d(n,i,j);
66          }
67        }
68    }
69    else
70    {
71      boundsLonSrc = domainSrc_->bounds_lon_1d;
72      boundsLatSrc = domainSrc_->bounds_lat_1d;
73    }
74  }
75  else // if domain source is rectilinear, not do anything now
76  {
77    nVertexSrc = constNVertex;
78  }
79
80  int localDomainDestSize = domainDest_->i_index.numElements();
81  int niDest = domainDest_->ni.getValue(), njDest = domainDest_->nj.getValue();
82  bool hasBoundDest = domainDest_->hasBounds;
83  if (hasBoundDest) nVertexDest = domainDest_->nvertex.getValue();
84  CArray<double,2> boundsLonDest(nVertexDest,localDomainDestSize);
85  CArray<double,2> boundsLatDest(nVertexDest,localDomainDestSize);
86
87  if (CDomain::type_attr::rectilinear == domainDest_->type) dstPole[2] = 1;
88  if (hasBoundDest)
89  {
90    if (!domainDest_->bounds_lon_2d.isEmpty())
91    {
92      for (j = 0; j < njDest; ++j)
93        for (i = 0; i < niDest; ++i)
94        {
95          k=j*niDest+i;
96          for(int n=0;n<nVertexDest;++n)
97          {
98            boundsLonDest(n,k) = domainDest_->bounds_lon_2d(n,i,j);
99            boundsLatDest(n,k) = domainDest_->bounds_lat_2d(n,i,j);
100          }
101        }
102    }
103    else
104    {
105      boundsLonDest = domainDest_->bounds_lon_1d;
106      boundsLatDest = domainDest_->bounds_lat_1d;
107    }
108  }
109  else
110  {
111    // Ok, fill in boundary values for rectangular domain
112    domainDest_->fillInRectilinearBoundLonLat(boundsLonDest, boundsLatDest);
113    nVertexDest = constNVertex;
114  }
115
116  // Ok, now use mapper to calculate
117  int nSrcLocal = domainSrc_->i_index.numElements();
118  int nDstLocal = domainDest_->i_index.numElements();
119  Mapper mapper(client->intraComm);
120  mapper.setVerbosity(PROGRESS) ;
121  mapper.setSourceMesh(boundsLonSrc.dataFirst(), boundsLatSrc.dataFirst(), nVertexSrc, nSrcLocal, &srcPole[0]);
122  mapper.setTargetMesh(boundsLonDest.dataFirst(), boundsLatDest.dataFirst(), nVertexDest, nDstLocal, &dstPole[0]);
123  std::vector<double> timings = mapper.computeWeights(orderInterp);
124
125  for (int idx = 0;  idx < mapper.nWeights; ++idx)
126  {
127    transformationMapping_[mapper.targetWeightId[idx]].push_back(mapper.sourceWeightId[idx]);
128    transformationWeight_[mapper.targetWeightId[idx]].push_back(mapper.remapMatrix[idx]);
129  }
130}
131
132/*!
133  Compute the index mapping between domain on grid source and one on grid destination
134*/
135void CDomainAlgorithmInterpolate::computeIndexSourceMapping()
136{
137  if (!interpDomain_->file.isEmpty())
138    readRemapInfo();
139  else
140    computeRemap();
141}
142
143/*!
144  Read remap information from file then distribute it among clients
145*/
146void CDomainAlgorithmInterpolate::readRemapInfo()
147{
148  CContext* context = CContext::getCurrent();
149  CContextClient* client=context->client;
150  int clientRank = client->clientRank;
151
152  std::string filename = interpDomain_->file.getValue();
153  std::map<int,std::vector<std::pair<int,double> > > interpMapValue;
154  readInterpolationInfo(filename, interpMapValue);
155
156  boost::unordered_map<size_t,int> globalIndexOfDomainDest;
157  int ni = domainDest_->ni.getValue();
158  int nj = domainDest_->nj.getValue();
159  int ni_glo = domainDest_->ni_glo.getValue();
160  size_t globalIndex;
161  int nIndexSize = domainDest_->i_index.numElements(), i_ind, j_ind;
162  for (int idx = 0; idx < nIndexSize; ++idx)
163  {
164    i_ind=domainDest_->i_index(idx) ;
165    j_ind=domainDest_->j_index(idx) ;
166
167    globalIndex = i_ind + j_ind * ni_glo;
168    globalIndexOfDomainDest[globalIndex] = clientRank;
169  }
170
171  CClientServerMappingDistributed domainIndexClientClientMapping(globalIndexOfDomainDest,
172                                                                 client->intraComm,
173                                                                 true);
174  CArray<size_t,1> globalIndexInterp(interpMapValue.size());
175  std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it,
176                                                                     ite = interpMapValue.end();
177  size_t globalIndexCount = 0;
178  for (it = itb; it != ite; ++it)
179  {
180    globalIndexInterp(globalIndexCount) = it->first;
181    ++globalIndexCount;
182  }
183
184  domainIndexClientClientMapping.computeServerIndexMapping(globalIndexInterp);
185  const std::map<int, std::vector<size_t> >& globalIndexInterpSendToClient = domainIndexClientClientMapping.getGlobalIndexOnServer();
186
187  //Inform each client number of index they will receive
188  int nbClient = client->clientSize;
189  int* sendBuff = new int[nbClient];
190  int* recvBuff = new int[nbClient];
191  for (int i = 0; i < nbClient; ++i) sendBuff[i] = 0;
192  int sendBuffSize = 0;
193  std::map<int, std::vector<size_t> >::const_iterator itbMap = globalIndexInterpSendToClient.begin(), itMap,
194                                                      iteMap = globalIndexInterpSendToClient.end();
195  for (itMap = itbMap; itMap != iteMap; ++itMap)
196  {
197    int sizeIndex = 0, mapSize = (itMap->second).size();
198    for (int idx = 0; idx < mapSize; ++idx)
199    {
200      sizeIndex += interpMapValue[(itMap->second)[idx]].size();
201    }
202    sendBuff[itMap->first] = sizeIndex;
203    sendBuffSize += sizeIndex;
204  }
205
206  MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm);
207
208  int* sendIndexDestBuff = new int [sendBuffSize];
209  int* sendIndexSrcBuff  = new int [sendBuffSize];
210  double* sendWeightBuff = new double [sendBuffSize];
211
212  std::vector<MPI_Request> sendRequest;
213  // Now send index and weight
214  int sendOffSet = 0;
215  for (itMap = itbMap; itMap != iteMap; ++itMap)
216  {
217    int k = 0;
218    int mapSize = (itMap->second).size();
219
220    for (int idx = 0; idx < mapSize; ++idx)
221    {
222      std::vector<std::pair<int,double> >& interpMap = interpMapValue[(itMap->second)[idx]];
223      for (int i = 0; i < interpMap.size(); ++i)
224      {
225        sendIndexDestBuff[k] = (itMap->second)[idx];
226        sendIndexSrcBuff[k]  = interpMap[i].first;
227        sendWeightBuff[k]    = interpMap[i].second;
228        ++k;
229      }
230    }
231
232    sendRequest.push_back(MPI_Request());
233    MPI_Isend(sendIndexDestBuff + sendOffSet,
234             k,
235             MPI_INT,
236             itMap->first,
237             7,
238             client->intraComm,
239             &sendRequest.back());
240    sendRequest.push_back(MPI_Request());
241    MPI_Isend(sendIndexSrcBuff + sendOffSet,
242             k,
243             MPI_INT,
244             itMap->first,
245             8,
246             client->intraComm,
247             &sendRequest.back());
248    sendRequest.push_back(MPI_Request());
249    MPI_Isend(sendWeightBuff + sendOffSet,
250             k,
251             MPI_DOUBLE,
252             itMap->first,
253             9,
254             client->intraComm,
255             &sendRequest.back());
256    sendOffSet += k;
257  }
258
259  int recvBuffSize = recvBuff[clientRank];
260  int* recvIndexDestBuff = new int [recvBuffSize];
261  int* recvIndexSrcBuff  = new int [recvBuffSize];
262  double* recvWeightBuff = new double [recvBuffSize];
263  int receivedSize = 0;
264  int clientSrcRank;
265  while (receivedSize < recvBuffSize)
266  {
267    MPI_Status recvStatus;
268    MPI_Recv((recvIndexDestBuff + receivedSize),
269             recvBuffSize,
270             MPI_INT,
271             MPI_ANY_SOURCE,
272             7,
273             client->intraComm,
274             &recvStatus);
275
276    int countBuff = 0;
277    MPI_Get_count(&recvStatus, MPI_INT, &countBuff);
278    clientSrcRank = recvStatus.MPI_SOURCE;
279
280    MPI_Recv((recvIndexSrcBuff + receivedSize),
281             recvBuffSize,
282             MPI_INT,
283             clientSrcRank,
284             8,
285             client->intraComm,
286             &recvStatus);
287
288    MPI_Recv((recvWeightBuff + receivedSize),
289             recvBuffSize,
290             MPI_DOUBLE,
291             clientSrcRank,
292             9,
293             client->intraComm,
294             &recvStatus);
295
296    for (int idx = 0; idx < countBuff; ++idx)
297    {
298      transformationMapping_[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvIndexSrcBuff + receivedSize + idx));
299      transformationWeight_[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvWeightBuff + receivedSize + idx));
300    }
301    receivedSize += countBuff;
302  }
303
304  std::vector<MPI_Status> requestStatus(sendRequest.size());
305  MPI_Wait(&sendRequest[0], &requestStatus[0]);
306
307  delete [] sendIndexDestBuff;
308  delete [] sendIndexSrcBuff;
309  delete [] sendWeightBuff;
310  delete [] recvIndexDestBuff;
311  delete [] recvIndexSrcBuff;
312  delete [] recvWeightBuff;
313  delete [] sendBuff;
314  delete [] recvBuff;
315}
316
317/*!
318  Read interpolation information from a file
319  \param [in] filename interpolation file
320  \param [in/out] interpMapValue Mapping between (global) index of domain on grid destination and
321         corresponding global index of domain and associated weight value on grid source
322*/
323void CDomainAlgorithmInterpolate::readInterpolationInfo(std::string& filename,
324                                                                std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
325{
326  int ncid ;
327  int weightDimId ;
328  size_t nbWeightGlo ;
329
330  CContext* context = CContext::getCurrent();
331  CContextClient* client=context->client;
332  int clientRank = client->clientRank;
333  int clientSize = client->clientSize;
334
335  nc_open(filename.c_str(),NC_NOWRITE, &ncid) ;
336  nc_inq_dimid(ncid,"n_weight",&weightDimId) ;
337  nc_inq_dimlen(ncid,weightDimId,&nbWeightGlo) ;
338
339  size_t nbWeight ;
340  size_t start ;
341  size_t div = nbWeightGlo/clientSize ;
342  size_t mod = nbWeightGlo%clientSize ;
343  if (clientRank < mod)
344  {
345    nbWeight=div+1 ;
346    start=clientRank*(div+1) ;
347  }
348  else
349  {
350    nbWeight=div ;
351    start= mod * (div+1) + (clientRank-mod) * div ;
352  }
353
354  double* weight=new double[nbWeight] ;
355  int weightId ;
356  nc_inq_varid (ncid, "weight", &weightId) ;
357  nc_get_vara_double(ncid, weightId, &start, &nbWeight, weight) ;
358
359  long* srcIndex=new long[nbWeight] ;
360  int srcIndexId ;
361  nc_inq_varid (ncid, "src_idx", &srcIndexId) ;
362  nc_get_vara_long(ncid, srcIndexId, &start, &nbWeight, srcIndex) ;
363
364  long* dstIndex=new long[nbWeight] ;
365  int dstIndexId ;
366  nc_inq_varid (ncid, "dst_idx", &dstIndexId) ;
367  nc_get_vara_long(ncid, dstIndexId, &start, &nbWeight, dstIndex) ;
368
369  for(size_t ind=0; ind<nbWeight;++ind)
370    interpMapValue[dstIndex[ind]-1].push_back(make_pair(srcIndex[ind]-1,weight[ind]));
371}
372
373}
Note: See TracBrowser for help on using the repository browser.