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

Last change on this file since 915 was 915, checked in by rlacroix, 8 years ago

Cosmetics: Fix intentation.

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