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

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

Implementing the reading of attributes of an axis from a file

+) 3d grid can be read directly from a file
+) Clean some redundant codes
+) Add new attribute declaration that allows to output only desired attributes

Test
+) On Curie
+) test_remap passes and result is correct

File size: 18.8 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  const double poleValue = 90.0;
42  const int constNVertex = 4; // Value by default number of vertex for rectangular domain
43  int nVertexSrc, nVertexDest;
44  nVertexSrc = nVertexDest = constNVertex;
45
46  // First of all, try to retrieve the boundary values of domain source and domain destination
47  int localDomainSrcSize = domainSrc_->i_index.numElements();
48  int niSrc = domainSrc_->ni.getValue(), njSrc = domainSrc_->nj.getValue();
49  bool hasBoundSrc = domainSrc_->hasBounds;
50  if (hasBoundSrc) nVertexSrc = domainSrc_->nvertex.getValue();
51  CArray<double,2> boundsLonSrc(nVertexSrc,localDomainSrcSize);
52  CArray<double,2> boundsLatSrc(nVertexSrc,localDomainSrcSize);
53
54  if (CDomain::type_attr::rectilinear == domainSrc_->type) srcPole[2] = 1;
55  if (hasBoundSrc)  // Suppose that domain source is curvilinear or unstructured
56  {
57    if (!domainSrc_->bounds_lon_2d.isEmpty())
58    {
59      for (j = 0; j < njSrc; ++j)
60        for (i = 0; i < niSrc; ++i)
61        {
62          k=j*niSrc+i;
63          for(int n=0;n<nVertexSrc;++n)
64          {
65            boundsLonSrc(n,k) = domainSrc_->bounds_lon_2d(n,i,j);
66            boundsLatSrc(n,k) = domainSrc_->bounds_lat_2d(n,i,j);
67          }
68        }
69    }
70    else
71    {
72      boundsLonSrc = domainSrc_->bounds_lon_1d;
73      boundsLatSrc = domainSrc_->bounds_lat_1d;
74    }
75  }
76  else // if domain source is rectilinear, not do anything now
77  {
78    bool isNorthPole = false;
79    bool isSouthPole = false;
80    if (std::abs(poleValue - std::abs(domainSrc_->lat_start)) < NumTraits<double>::epsilon()) isNorthPole = true;
81    if (std::abs(poleValue - std::abs(domainSrc_->lat_end)) < NumTraits<double>::epsilon()) isSouthPole = true;
82
83    nVertexSrc = constNVertex;
84    domainSrc_->fillInRectilinearBoundLonLat(boundsLonSrc, boundsLatSrc, isNorthPole, isSouthPole);
85  }
86
87  std::map<int,std::vector<std::pair<int,double> > > interpMapValueNorthPole;
88  std::map<int,std::vector<std::pair<int,double> > > interpMapValueSouthPole;
89
90  int localDomainDestSize = domainDest_->i_index.numElements();
91  int niDest = domainDest_->ni.getValue(), njDest = domainDest_->nj.getValue();
92  bool hasBoundDest = domainDest_->hasBounds;
93  if (hasBoundDest) nVertexDest = domainDest_->nvertex.getValue();
94  CArray<double,2> boundsLonDest(nVertexDest,localDomainDestSize);
95  CArray<double,2> boundsLatDest(nVertexDest,localDomainDestSize);
96
97  if (CDomain::type_attr::rectilinear == domainDest_->type) dstPole[2] = 1;
98  if (hasBoundDest)
99  {
100    if (!domainDest_->bounds_lon_2d.isEmpty())
101    {
102      for (j = 0; j < njDest; ++j)
103        for (i = 0; i < niDest; ++i)
104        {
105          k=j*niDest+i;
106          for(int n=0;n<nVertexDest;++n)
107          {
108            boundsLonDest(n,k) = domainDest_->bounds_lon_2d(n,i,j);
109            boundsLatDest(n,k) = domainDest_->bounds_lat_2d(n,i,j);
110          }
111        }
112    }
113    else
114    {
115      boundsLonDest = domainDest_->bounds_lon_1d;
116      boundsLatDest = domainDest_->bounds_lat_1d;
117    }
118  }
119  else
120  {
121    bool isNorthPole = false;
122    bool isSouthPole = false;
123    if (std::abs(poleValue - std::abs(domainDest_->lat_start)) < NumTraits<double>::epsilon()) isNorthPole = true;
124    if (std::abs(poleValue - std::abs(domainDest_->lat_end)) < NumTraits<double>::epsilon()) isSouthPole = true;
125    if (isNorthPole && (0 == domainDest_->jbegin.getValue()))
126    {
127      int ibegin = domainDest_->ibegin.getValue();
128      for (i = 0; i < niDest; ++i)
129      {
130        interpMapValueNorthPole[i+ibegin];
131      }
132    }
133
134    if (isSouthPole && (domainDest_->nj_glo.getValue() == (domainDest_->jbegin.getValue() + njDest)))
135    {
136      int ibegin = domainDest_->ibegin.getValue();
137      int njGlo = domainDest_->nj_glo.getValue();
138      int niGlo = domainDest_->ni_glo.getValue();
139      for (i = 0; i < niDest; ++i)
140      {
141        k = (njGlo - 1)*niGlo + i + ibegin;
142        interpMapValueSouthPole[k];
143      }
144    }
145
146    // Ok, fill in boundary values for rectangular domain
147    nVertexDest = constNVertex;
148    domainDest_->fillInRectilinearBoundLonLat(boundsLonDest, boundsLatDest, isNorthPole, isSouthPole);
149  }
150
151
152
153  // Ok, now use mapper to calculate
154  int nSrcLocal = domainSrc_->i_index.numElements();
155  int nDstLocal = domainDest_->i_index.numElements();
156  long int * globalSrc = new long int [nSrcLocal];
157  long int * globalDst = new long int [nDstLocal];
158
159  long int globalIndex;
160  int i_ind, j_ind;
161  for (int idx = 0; idx < nSrcLocal; ++idx)
162  {
163    i_ind=domainSrc_->i_index(idx) ;
164    j_ind=domainSrc_->j_index(idx) ;
165
166    globalIndex = i_ind + j_ind * domainSrc_->ni_glo;
167    globalSrc[idx] = globalIndex;
168  }
169
170  for (int idx = 0; idx < nDstLocal; ++idx)
171  {
172    i_ind=domainDest_->i_index(idx) ;
173    j_ind=domainDest_->j_index(idx) ;
174
175    globalIndex = i_ind + j_ind * domainDest_->ni_glo;
176    globalDst[idx] = globalIndex;
177  }
178
179
180  // Calculate weight index
181  Mapper mapper(client->intraComm);
182  mapper.setVerbosity(PROGRESS) ;
183  mapper.setSourceMesh(boundsLonSrc.dataFirst(), boundsLatSrc.dataFirst(), nVertexSrc, nSrcLocal, &srcPole[0], globalSrc);
184  mapper.setTargetMesh(boundsLonDest.dataFirst(), boundsLatDest.dataFirst(), nVertexDest, nDstLocal, &dstPole[0], globalDst);
185  std::vector<double> timings = mapper.computeWeights(orderInterp);
186
187  std::map<int,std::vector<std::pair<int,double> > > interpMapValue;
188  std::map<int,std::vector<std::pair<int,double> > >::const_iterator iteNorthPole = interpMapValueNorthPole.end(),
189                                                                     iteSouthPole = interpMapValueSouthPole.end();
190  for (int idx = 0;  idx < mapper.nWeights; ++idx)
191  {
192    interpMapValue[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
193    if (iteNorthPole != interpMapValueNorthPole.find(mapper.targetWeightId[idx]))
194    {
195      interpMapValueNorthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
196    }
197
198    if (iteSouthPole != interpMapValueSouthPole.find(mapper.targetWeightId[idx]))
199    {
200      interpMapValueSouthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
201    }
202  }
203  int niGloDst = domainDest_->ni_glo.getValue();
204  processPole(interpMapValueNorthPole, niGloDst);
205  processPole(interpMapValueSouthPole, niGloDst);
206
207  if (!interpMapValueNorthPole.empty())
208  {
209     std::map<int,std::vector<std::pair<int,double> > >::iterator itNorthPole = interpMapValueNorthPole.begin();
210     for (; itNorthPole != iteNorthPole; ++itNorthPole)
211     {
212       if (!(itNorthPole->second.empty()))
213        itNorthPole->second.swap(interpMapValue[itNorthPole->first]);
214     }
215  }
216
217  if (!interpMapValueSouthPole.empty())
218  {
219     std::map<int,std::vector<std::pair<int,double> > >::iterator itSouthPole = interpMapValueSouthPole.begin();
220     for (; itSouthPole != iteSouthPole; ++itSouthPole)
221     {
222       if (!(itSouthPole->second.empty()))
223        itSouthPole->second.swap(interpMapValue[itSouthPole->first]);
224     }
225  }
226
227  exchangeRemapInfo(interpMapValue);
228
229  delete [] globalSrc;
230  delete [] globalDst;
231}
232
233void CDomainAlgorithmInterpolate::processPole(std::map<int,std::vector<std::pair<int,double> > >& interMapValuePole,
234                                              int nbGlobalPointOnPole)
235{
236  CContext* context = CContext::getCurrent();
237  CContextClient* client=context->client;
238
239  MPI_Comm poleComme(MPI_COMM_NULL);
240  MPI_Comm_split(client->intraComm, interMapValuePole.empty() ? MPI_UNDEFINED : 1, 0, &poleComme);
241  if (MPI_COMM_NULL != poleComme)
242  {
243    int nbClientPole;
244    MPI_Comm_size(poleComme, &nbClientPole);
245
246    std::map<int,std::vector<std::pair<int,double> > >::iterator itePole = interMapValuePole.end(), itPole,
247                                                                 itbPole = interMapValuePole.begin();
248
249    int nbWeight = 0;
250    for (itPole = itbPole; itPole != itePole; ++itPole)
251       nbWeight += itPole->second.size();
252
253    std::vector<int> recvCount(nbClientPole,0);
254    std::vector<int> displ(nbClientPole,0);
255    MPI_Allgather(&nbWeight,1,MPI_INT,&recvCount[0],1,MPI_INT,poleComme) ;
256
257    displ[0]=0;
258    for(int n=1;n<nbClientPole;++n) displ[n]=displ[n-1]+recvCount[n-1] ;
259    int recvSize=displ[nbClientPole-1]+recvCount[nbClientPole-1] ;
260
261    std::vector<int> sendSourceIndexBuff(nbWeight);
262    std::vector<double> sendSourceWeightBuff(nbWeight);
263    int k = 0;
264    for (itPole = itbPole; itPole != itePole; ++itPole)
265    {
266      for (int idx = 0; idx < itPole->second.size(); ++idx)
267      {
268        sendSourceIndexBuff[k] = (itPole->second)[idx].first;
269        sendSourceWeightBuff[k] = (itPole->second)[idx].second;
270        ++k;
271      }
272    }
273
274    std::vector<int> recvSourceIndexBuff(recvSize);
275    std::vector<double> recvSourceWeightBuff(recvSize);
276
277    // Gather all index and weight for pole
278    MPI_Allgatherv(&sendSourceIndexBuff[0],nbWeight,MPI_INT,&recvSourceIndexBuff[0],&recvCount[0],&displ[0],MPI_INT,poleComme);
279    MPI_Allgatherv(&sendSourceWeightBuff[0],nbWeight,MPI_DOUBLE,&recvSourceWeightBuff[0],&recvCount[0],&displ[0],MPI_DOUBLE,poleComme);
280
281    std::map<int,double> recvTemp;
282    for (int idx = 0; idx < recvSize; ++idx)
283    {
284      if (recvTemp.end() != recvTemp.find(recvSourceIndexBuff[idx]))
285        recvTemp[recvSourceIndexBuff[idx]] += recvSourceWeightBuff[idx]/nbGlobalPointOnPole;
286      else
287        recvTemp[recvSourceIndexBuff[idx]] = 0.0;
288    }
289
290    std::map<int,double>::const_iterator itRecvTemp, itbRecvTemp = recvTemp.begin(), iteRecvTemp = recvTemp.end();
291
292    for (itPole = itbPole; itPole != itePole; ++itPole)
293    {
294      itPole->second.clear();
295      for (itRecvTemp = itbRecvTemp; itRecvTemp != iteRecvTemp; ++itRecvTemp)
296          itPole->second.push_back(make_pair(itRecvTemp->first, itRecvTemp->second));
297    }
298  }
299
300}
301
302/*!
303  Compute the index mapping between domain on grid source and one on grid destination
304*/
305void CDomainAlgorithmInterpolate::computeIndexSourceMapping()
306{
307  if (!interpDomain_->file.isEmpty())
308    readRemapInfo();
309  else
310    computeRemap();
311}
312
313void CDomainAlgorithmInterpolate::readRemapInfo()
314{
315  CContext* context = CContext::getCurrent();
316  CContextClient* client=context->client;
317  int clientRank = client->clientRank;
318
319  std::string filename = interpDomain_->file.getValue();
320  std::map<int,std::vector<std::pair<int,double> > > interpMapValue;
321  readInterpolationInfo(filename, interpMapValue);
322
323  exchangeRemapInfo(interpMapValue);
324}
325
326
327/*!
328  Read remap information from file then distribute it among clients
329*/
330void CDomainAlgorithmInterpolate::exchangeRemapInfo(const std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
331{
332  CContext* context = CContext::getCurrent();
333  CContextClient* client=context->client;
334  int clientRank = client->clientRank;
335
336  boost::unordered_map<size_t,int> globalIndexOfDomainDest;
337  int ni = domainDest_->ni.getValue();
338  int nj = domainDest_->nj.getValue();
339  int ni_glo = domainDest_->ni_glo.getValue();
340  size_t globalIndex;
341  int nIndexSize = domainDest_->i_index.numElements(), i_ind, j_ind;
342  for (int idx = 0; idx < nIndexSize; ++idx)
343  {
344    i_ind=domainDest_->i_index(idx) ;
345    j_ind=domainDest_->j_index(idx) ;
346
347    globalIndex = i_ind + j_ind * ni_glo;
348    globalIndexOfDomainDest[globalIndex] = clientRank;
349  }
350
351  CClientServerMappingDistributed domainIndexClientClientMapping(globalIndexOfDomainDest,
352                                                                 client->intraComm,
353                                                                 true);
354  CArray<size_t,1> globalIndexInterp(interpMapValue.size());
355  std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it,
356                                                                     ite = interpMapValue.end();
357  size_t globalIndexCount = 0;
358  for (it = itb; it != ite; ++it)
359  {
360    globalIndexInterp(globalIndexCount) = it->first;
361    ++globalIndexCount;
362  }
363
364  domainIndexClientClientMapping.computeServerIndexMapping(globalIndexInterp);
365  const std::map<int, std::vector<size_t> >& globalIndexInterpSendToClient = domainIndexClientClientMapping.getGlobalIndexOnServer();
366
367  //Inform each client number of index they will receive
368  int nbClient = client->clientSize;
369  int* sendBuff = new int[nbClient];
370  int* recvBuff = new int[nbClient];
371  for (int i = 0; i < nbClient; ++i)
372  {
373    sendBuff[i] = 0;
374    recvBuff[i] = 0;
375  }
376  int sendBuffSize = 0;
377  std::map<int, std::vector<size_t> >::const_iterator itbMap = globalIndexInterpSendToClient.begin(), itMap,
378                                                      iteMap = globalIndexInterpSendToClient.end();
379  for (itMap = itbMap; itMap != iteMap; ++itMap)
380  {
381    const std::vector<size_t>& tmp = itMap->second;
382    int sizeIndex = 0, mapSize = (itMap->second).size();
383    for (int idx = 0; idx < mapSize; ++idx)
384    {
385      sizeIndex += interpMapValue.at((itMap->second)[idx]).size();
386    }
387    sendBuff[itMap->first] = sizeIndex;
388    sendBuffSize += sizeIndex;
389  }
390
391
392  MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm);
393
394  int* sendIndexDestBuff = new int [sendBuffSize];
395  int* sendIndexSrcBuff  = new int [sendBuffSize];
396  double* sendWeightBuff = new double [sendBuffSize];
397
398  std::vector<MPI_Request> sendRequest;
399
400  int sendOffSet = 0, l = 0;
401  for (itMap = itbMap; itMap != iteMap; ++itMap)
402  {
403    const std::vector<size_t>& indexToSend = itMap->second;
404    int mapSize = indexToSend.size();
405    int k = 0;
406    for (int idx = 0; idx < mapSize; ++idx)
407    {
408      const std::vector<std::pair<int,double> >& interpMap = interpMapValue.at(indexToSend[idx]);
409      for (int i = 0; i < interpMap.size(); ++i)
410      {
411        sendIndexDestBuff[l] = indexToSend[idx];
412        sendIndexSrcBuff[l]  = interpMap[i].first;
413        sendWeightBuff[l]    = interpMap[i].second;
414        ++k;
415        ++l;
416      }
417    }
418
419    sendRequest.push_back(MPI_Request());
420    MPI_Isend(sendIndexDestBuff + sendOffSet,
421             k,
422             MPI_INT,
423             itMap->first,
424             7,
425             client->intraComm,
426             &sendRequest.back());
427    sendRequest.push_back(MPI_Request());
428    MPI_Isend(sendIndexSrcBuff + sendOffSet,
429             k,
430             MPI_INT,
431             itMap->first,
432             8,
433             client->intraComm,
434             &sendRequest.back());
435    sendRequest.push_back(MPI_Request());
436    MPI_Isend(sendWeightBuff + sendOffSet,
437             k,
438             MPI_DOUBLE,
439             itMap->first,
440             9,
441             client->intraComm,
442             &sendRequest.back());
443    sendOffSet += k;
444  }
445
446  int recvBuffSize = recvBuff[clientRank];
447  int* recvIndexDestBuff = new int [recvBuffSize];
448  int* recvIndexSrcBuff  = new int [recvBuffSize];
449  double* recvWeightBuff = new double [recvBuffSize];
450  int receivedSize = 0;
451  int clientSrcRank;
452  while (receivedSize < recvBuffSize)
453  {
454    MPI_Status recvStatus;
455    MPI_Recv((recvIndexDestBuff + receivedSize),
456             recvBuffSize,
457             MPI_INT,
458             MPI_ANY_SOURCE,
459             7,
460             client->intraComm,
461             &recvStatus);
462
463    int countBuff = 0;
464    MPI_Get_count(&recvStatus, MPI_INT, &countBuff);
465    clientSrcRank = recvStatus.MPI_SOURCE;
466
467    MPI_Recv((recvIndexSrcBuff + receivedSize),
468             recvBuffSize,
469             MPI_INT,
470             clientSrcRank,
471             8,
472             client->intraComm,
473             &recvStatus);
474
475    MPI_Recv((recvWeightBuff + receivedSize),
476             recvBuffSize,
477             MPI_DOUBLE,
478             clientSrcRank,
479             9,
480             client->intraComm,
481             &recvStatus);
482
483    for (int idx = 0; idx < countBuff; ++idx)
484    {
485      transformationMapping_[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvIndexSrcBuff + receivedSize + idx));
486      transformationWeight_[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvWeightBuff + receivedSize + idx));
487    }
488    receivedSize += countBuff;
489  }
490
491  std::vector<MPI_Status> requestStatus(sendRequest.size());
492  MPI_Waitall(sendRequest.size(), &sendRequest[0], MPI_STATUS_IGNORE);
493
494  delete [] sendIndexDestBuff;
495  delete [] sendIndexSrcBuff;
496  delete [] sendWeightBuff;
497  delete [] recvIndexDestBuff;
498  delete [] recvIndexSrcBuff;
499  delete [] recvWeightBuff;
500  delete [] sendBuff;
501  delete [] recvBuff;
502}
503
504/*!
505  Read interpolation information from a file
506  \param [in] filename interpolation file
507  \param [in/out] interpMapValue Mapping between (global) index of domain on grid destination and
508         corresponding global index of domain and associated weight value on grid source
509*/
510void CDomainAlgorithmInterpolate::readInterpolationInfo(std::string& filename,
511                                                        std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
512{
513  int ncid ;
514  int weightDimId ;
515  size_t nbWeightGlo ;
516
517  CContext* context = CContext::getCurrent();
518  CContextClient* client=context->client;
519  int clientRank = client->clientRank;
520  int clientSize = client->clientSize;
521
522  nc_open(filename.c_str(),NC_NOWRITE, &ncid) ;
523  nc_inq_dimid(ncid,"n_weight",&weightDimId) ;
524  nc_inq_dimlen(ncid,weightDimId,&nbWeightGlo) ;
525
526  size_t nbWeight ;
527  size_t start ;
528  size_t div = nbWeightGlo/clientSize ;
529  size_t mod = nbWeightGlo%clientSize ;
530  if (clientRank < mod)
531  {
532    nbWeight=div+1 ;
533    start=clientRank*(div+1) ;
534  }
535  else
536  {
537    nbWeight=div ;
538    start= mod * (div+1) + (clientRank-mod) * div ;
539  }
540
541  double* weight=new double[nbWeight] ;
542  int weightId ;
543  nc_inq_varid (ncid, "weight", &weightId) ;
544  nc_get_vara_double(ncid, weightId, &start, &nbWeight, weight) ;
545
546  long* srcIndex=new long[nbWeight] ;
547  int srcIndexId ;
548  nc_inq_varid (ncid, "src_idx", &srcIndexId) ;
549  nc_get_vara_long(ncid, srcIndexId, &start, &nbWeight, srcIndex) ;
550
551  long* dstIndex=new long[nbWeight] ;
552  int dstIndexId ;
553  nc_inq_varid (ncid, "dst_idx", &dstIndexId) ;
554  nc_get_vara_long(ncid, dstIndexId, &start, &nbWeight, dstIndex) ;
555
556  for(size_t ind=0; ind<nbWeight;++ind)
557    interpMapValue[dstIndex[ind]-1].push_back(make_pair(srcIndex[ind]-1,weight[ind]));
558}
559
560}
Note: See TracBrowser for help on using the repository browser.