source: XIOS/dev/branch_openmp/src/transformation/domain_algorithm_interpolate.cpp @ 1328

Last change on this file since 1328 was 1328, checked in by yushan, 6 years ago

dev_omp

File size: 31.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#include "mpi_tag.hpp"
18#include "domain.hpp"
19#include "grid_transformation_factory_impl.hpp"
20#include "interpolate_domain.hpp"
21#include "grid.hpp"
22using namespace ep_lib;
23
24namespace xios {
25CGenericAlgorithmTransformation* CDomainAlgorithmInterpolate::create(CGrid* gridDst, CGrid* gridSrc,
26                                                                     CTransformation<CDomain>* transformation,
27                                                                     int elementPositionInGrid,
28                                                                     std::map<int, int>& elementPositionInGridSrc2ScalarPosition,
29                                                                     std::map<int, int>& elementPositionInGridSrc2AxisPosition,
30                                                                     std::map<int, int>& elementPositionInGridSrc2DomainPosition,
31                                                                     std::map<int, int>& elementPositionInGridDst2ScalarPosition,
32                                                                     std::map<int, int>& elementPositionInGridDst2AxisPosition,
33                                                                     std::map<int, int>& elementPositionInGridDst2DomainPosition)
34{
35  std::vector<CDomain*> domainListDestP = gridDst->getDomains();
36  std::vector<CDomain*> domainListSrcP  = gridSrc->getDomains();
37
38  CInterpolateDomain* interpolateDomain = dynamic_cast<CInterpolateDomain*> (transformation);
39  int domainDstIndex = elementPositionInGridDst2DomainPosition[elementPositionInGrid];
40  int domainSrcIndex = elementPositionInGridSrc2DomainPosition[elementPositionInGrid];
41
42  return (new CDomainAlgorithmInterpolate(domainListDestP[domainDstIndex], domainListSrcP[domainSrcIndex], interpolateDomain));
43}
44
45bool CDomainAlgorithmInterpolate::registerTrans()
46{
47  CGridTransformationFactory<CDomain>::registerTransformation(TRANS_INTERPOLATE_DOMAIN, create);
48}
49
50CDomainAlgorithmInterpolate::CDomainAlgorithmInterpolate(CDomain* domainDestination, CDomain* domainSource, CInterpolateDomain* interpDomain)
51: CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain), writeToFile_(false), readFromFile_(false)
52{
53  CContext* context = CContext::getCurrent();
54  interpDomain_->checkValid(domainSource);
55  fileToReadWrite_ = "xios_interpolation_weights_";
56
57  if (interpDomain_->weight_filename.isEmpty())
58  {
59    fileToReadWrite_ += context->getId() + "_" + 
60                    domainSource->getDomainOutputName() + "_" + 
61                    domainDestination->getDomainOutputName() + ".nc";   
62  }
63  else 
64    fileToReadWrite_ = interpDomain_->weight_filename;
65
66  ifstream f(fileToReadWrite_.c_str()); 
67  switch (interpDomain_->mode)
68  {
69    case CInterpolateDomain::mode_attr::read:
70      readFromFile_ = true;     
71      break;
72    case CInterpolateDomain::mode_attr::compute:
73      readFromFile_ = false;
74      break;
75    case CInterpolateDomain::mode_attr::read_or_compute:     
76      if (!f.good())
77        readFromFile_ = false;
78      else
79        readFromFile_ = true;
80      break;
81    default:
82      break;
83  } 
84
85  writeToFile_ = interpDomain_->write_weight; 
86   
87}
88
89/*!
90  Compute remap with integrated remap calculation module
91*/
92void CDomainAlgorithmInterpolate::computeRemap()
93{
94  using namespace sphereRemap;
95
96  CContext* context = CContext::getCurrent();
97  CContextClient* client=context->client;
98  int clientRank = client->clientRank;
99  int i, j, k, idx;
100  std::vector<double> srcPole(3,0), dstPole(3,0);
101  int orderInterp = interpDomain_->order.getValue();
102  bool renormalize ;
103  bool quantity ;
104
105  if (interpDomain_->renormalize.isEmpty()) renormalize=true;
106  else renormalize = interpDomain_->renormalize;
107
108  if (interpDomain_->quantity.isEmpty()) quantity=false;
109  else quantity = interpDomain_->quantity;
110
111  const double poleValue = 90.0;
112  const int constNVertex = 4; // Value by default number of vertex for rectangular domain
113  int nVertexSrc, nVertexDest;
114  nVertexSrc = nVertexDest = constNVertex;
115
116
117
118  // First of all, try to retrieve the boundary values of domain source and domain destination
119  int localDomainSrcSize = domainSrc_->i_index.numElements();
120  int niSrc = domainSrc_->ni.getValue(), njSrc = domainSrc_->nj.getValue();
121  bool hasBoundSrc = domainSrc_->hasBounds;
122  if (hasBoundSrc) nVertexSrc = domainSrc_->nvertex.getValue();
123  CArray<double,2> boundsLonSrc(nVertexSrc,localDomainSrcSize);
124  CArray<double,2> boundsLatSrc(nVertexSrc,localDomainSrcSize);
125
126  if (domainSrc_->hasPole) srcPole[2] = 1;
127  if (hasBoundSrc)  // Suppose that domain source is curvilinear or unstructured
128  {
129    if (!domainSrc_->bounds_lon_2d.isEmpty())
130    {
131      for (j = 0; j < njSrc; ++j)
132        for (i = 0; i < niSrc; ++i)
133        {
134          k=j*niSrc+i;
135          for(int n=0;n<nVertexSrc;++n)
136          {
137            boundsLonSrc(n,k) = domainSrc_->bounds_lon_2d(n,i,j);
138            boundsLatSrc(n,k) = domainSrc_->bounds_lat_2d(n,i,j);
139          }
140        }
141    }
142    else
143    {
144      boundsLonSrc = domainSrc_->bounds_lon_1d;
145      boundsLatSrc = domainSrc_->bounds_lat_1d;
146    }
147  }
148  else // if domain source is rectilinear, not do anything now
149  {
150    CArray<double,1> lon_g ;
151    CArray<double,1> lat_g ;
152
153    if (!domainSrc_->lonvalue_1d.isEmpty() && !domainSrc_->latvalue_1d.isEmpty())
154    {
155      domainSrc_->AllgatherRectilinearLonLat(domainSrc_->lonvalue_1d,domainSrc_->latvalue_1d, lon_g,lat_g) ;
156    }
157    else if (! domainSrc_->latvalue_rectilinear_read_from_file.isEmpty() && ! domainSrc_->lonvalue_rectilinear_read_from_file.isEmpty() )
158    {
159      lat_g=domainSrc_->latvalue_rectilinear_read_from_file ;
160      lon_g=domainSrc_->lonvalue_rectilinear_read_from_file ;
161    }
162    else if (!domainSrc_->lon_start.isEmpty() && !domainSrc_->lon_end.isEmpty() &&
163             !domainSrc_->lat_start.isEmpty() && !domainSrc_->lat_end.isEmpty())
164    {
165      double step=(domainSrc_->lon_end-domainSrc_->lon_start)/domainSrc_->ni_glo ;
166      for (int i=0; i<domainSrc_->ni_glo; ++i) lon_g(i)=domainSrc_->lon_start+i*step ;
167      step=(domainSrc_->lat_end-domainSrc_->lat_start)/domainSrc_->nj_glo ;
168      for (int i=0; i<domainSrc_->ni_glo; ++i) lat_g(i)=domainSrc_->lat_start+i*step ;
169    }
170    else ERROR("void CDomainAlgorithmInterpolate::computeRemap()",<<"Cannot compute bounds for rectilinear domain") ;
171
172    nVertexSrc = constNVertex;
173    domainSrc_->fillInRectilinearBoundLonLat(lon_g,lat_g, boundsLonSrc, boundsLatSrc);
174  }
175
176  std::map<int,std::vector<std::pair<int,double> > > interpMapValueNorthPole;
177  std::map<int,std::vector<std::pair<int,double> > > interpMapValueSouthPole;
178
179  int localDomainDestSize = domainDest_->i_index.numElements();
180  int niDest = domainDest_->ni.getValue(), njDest = domainDest_->nj.getValue();
181  bool hasBoundDest = domainDest_->hasBounds;
182  if (hasBoundDest) nVertexDest = domainDest_->nvertex.getValue();
183  CArray<double,2> boundsLonDest(nVertexDest,localDomainDestSize);
184  CArray<double,2> boundsLatDest(nVertexDest,localDomainDestSize);
185
186  if (domainDest_->hasPole) dstPole[2] = 1;
187  if (hasBoundDest)
188  {
189    if (!domainDest_->bounds_lon_2d.isEmpty())
190    {
191      for (j = 0; j < njDest; ++j)
192        for (i = 0; i < niDest; ++i)
193        {
194          k=j*niDest+i;
195          for(int n=0;n<nVertexDest;++n)
196          {
197            boundsLonDest(n,k) = domainDest_->bounds_lon_2d(n,i,j);
198            boundsLatDest(n,k) = domainDest_->bounds_lat_2d(n,i,j);
199          }
200        }
201    }
202    else
203    {
204      boundsLonDest = domainDest_->bounds_lon_1d;
205      boundsLatDest = domainDest_->bounds_lat_1d;
206    }
207  }
208  else
209  {
210    bool isNorthPole = false;
211    bool isSouthPole = false;
212
213    CArray<double,1> lon_g ;
214    CArray<double,1> lat_g ;
215
216    if (!domainDest_->lonvalue_1d.isEmpty() && !domainDest_->latvalue_1d.isEmpty())
217    {
218      domainDest_->AllgatherRectilinearLonLat(domainDest_->lonvalue_1d,domainDest_->latvalue_1d, lon_g,lat_g) ;
219    }
220    else if (! domainDest_->latvalue_rectilinear_read_from_file.isEmpty() && ! domainDest_->lonvalue_rectilinear_read_from_file.isEmpty() )
221    {
222      lat_g=domainDest_->latvalue_rectilinear_read_from_file ;
223      lon_g=domainDest_->lonvalue_rectilinear_read_from_file ;
224    }
225    else if (!domainDest_->lon_start.isEmpty() && !domainDest_->lon_end.isEmpty() &&
226             !domainDest_->lat_start.isEmpty() && !domainDest_->lat_end.isEmpty())
227    {
228      double step=(domainDest_->lon_end-domainDest_->lon_start)/domainDest_->ni_glo ;
229      for(int i=0; i<domainDest_->ni_glo; ++i) lon_g(i)=domainDest_->lon_start+i*step ;
230      step=(domainDest_->lat_end-domainDest_->lat_start)/domainDest_->nj_glo ;
231      for(int i=0; i<domainDest_->ni_glo; ++i) lat_g(i)=domainDest_->lat_start+i*step ;
232    }
233    else ERROR("void CDomainAlgorithmInterpolate::computeRemap()",<<"Cannot compute bounds for rectilinear domain") ;
234   
235    if (std::abs(poleValue - std::abs(lat_g(0))) < NumTraits<double>::epsilon()) isNorthPole = true;
236    if (std::abs(poleValue - std::abs(lat_g(domainDest_->nj_glo-1))) < NumTraits<double>::epsilon()) isSouthPole = true;
237
238
239
240
241    if (isNorthPole && (0 == domainDest_->jbegin.getValue()))
242    {
243      int ibegin = domainDest_->ibegin.getValue();
244      for (i = 0; i < niDest; ++i)
245      {
246        interpMapValueNorthPole[i+ibegin];
247      }
248    }
249
250    if (isSouthPole && (domainDest_->nj_glo.getValue() == (domainDest_->jbegin.getValue() + njDest)))
251    {
252      int ibegin = domainDest_->ibegin.getValue();
253      int njGlo = domainDest_->nj_glo.getValue();
254      int niGlo = domainDest_->ni_glo.getValue();
255      for (i = 0; i < niDest; ++i)
256      {
257        k = (njGlo - 1)*niGlo + i + ibegin;
258        interpMapValueSouthPole[k];
259      }
260    }
261
262    // Ok, fill in boundary values for rectangular domain
263    nVertexDest = constNVertex;
264    domainDest_->fillInRectilinearBoundLonLat(lon_g,lat_g, boundsLonDest, boundsLatDest);
265  }
266
267
268
269  // Ok, now use mapper to calculate
270  int nSrcLocal = domainSrc_->i_index.numElements();
271  int nDstLocal = domainDest_->i_index.numElements();
272  long int * globalSrc = new long int [nSrcLocal];
273  long int * globalDst = new long int [nDstLocal];
274
275  long int globalIndex;
276  int i_ind, j_ind;
277  for (int idx = 0; idx < nSrcLocal; ++idx)
278  {
279    i_ind=domainSrc_->i_index(idx) ;
280    j_ind=domainSrc_->j_index(idx) ;
281
282    globalIndex = i_ind + j_ind * domainSrc_->ni_glo;
283    globalSrc[idx] = globalIndex;
284  }
285
286
287  for (int idx = 0; idx < nDstLocal; ++idx)
288  {
289    i_ind=domainDest_->i_index(idx) ;
290    j_ind=domainDest_->j_index(idx) ;
291
292    globalIndex = i_ind + j_ind * domainDest_->ni_glo;
293    globalDst[idx] = globalIndex;
294  }
295
296
297  // Calculate weight index
298  Mapper mapper(client->intraComm);
299  mapper.setVerbosity(PROGRESS) ;
300
301     
302  // supress masked data for the source
303  int nSrcLocalUnmasked = 0 ;
304  for (int idx=0 ; idx < nSrcLocal; idx++) if (domainSrc_->localMask(idx)) ++nSrcLocalUnmasked ;
305
306
307  CArray<double,2> boundsLonSrcUnmasked(nVertexSrc,nSrcLocalUnmasked);
308  CArray<double,2> boundsLatSrcUnmasked(nVertexSrc,nSrcLocalUnmasked);
309  long int * globalSrcUnmasked = new long int [nSrcLocalUnmasked];
310
311  nSrcLocalUnmasked=0 ;
312  for (int idx=0 ; idx < nSrcLocal; idx++)
313  {
314    if (domainSrc_->localMask(idx))
315    {
316      for(int n=0;n<nVertexSrc;++n)
317      {
318        boundsLonSrcUnmasked(n,nSrcLocalUnmasked) = boundsLonSrc(n,idx) ;
319        boundsLatSrcUnmasked(n,nSrcLocalUnmasked) = boundsLatSrc(n,idx) ;
320      }
321      globalSrcUnmasked[nSrcLocalUnmasked]=globalSrc[idx] ;
322      ++nSrcLocalUnmasked ;
323    }
324  }
325
326
327  int nDstLocalUnmasked = 0 ;
328  for (int idx=0 ; idx < nDstLocal; idx++) if (domainDest_->localMask(idx)) ++nDstLocalUnmasked ;
329
330  CArray<double,2> boundsLonDestUnmasked(nVertexDest,nDstLocalUnmasked);
331  CArray<double,2> boundsLatDestUnmasked(nVertexDest,nDstLocalUnmasked);
332  long int * globalDstUnmasked = new long int [nDstLocalUnmasked];
333
334  nDstLocalUnmasked=0 ;
335  for (int idx=0 ; idx < nDstLocal; idx++)
336  {
337    if (domainDest_->localMask(idx))
338    {
339      for(int n=0;n<nVertexDest;++n)
340      {
341        boundsLonDestUnmasked(n,nDstLocalUnmasked) = boundsLonDest(n,idx) ;
342        boundsLatDestUnmasked(n,nDstLocalUnmasked) = boundsLatDest(n,idx) ;
343      }
344      globalDstUnmasked[nDstLocalUnmasked]=globalDst[idx] ;
345      ++nDstLocalUnmasked ;
346    }
347  }
348
349  mapper.setSourceMesh(boundsLonSrcUnmasked.dataFirst(), boundsLatSrcUnmasked.dataFirst(), nVertexSrc, nSrcLocalUnmasked, &srcPole[0], globalSrcUnmasked);
350  mapper.setTargetMesh(boundsLonDestUnmasked.dataFirst(), boundsLatDestUnmasked.dataFirst(), nVertexDest, nDstLocalUnmasked, &dstPole[0], globalDstUnmasked);
351
352  std::vector<double> timings = mapper.computeWeights(orderInterp,renormalize,quantity);
353
354  std::map<int,std::vector<std::pair<int,double> > > interpMapValue;
355  std::map<int,std::vector<std::pair<int,double> > >::const_iterator iteNorthPole = interpMapValueNorthPole.end(),
356                                                                     iteSouthPole = interpMapValueSouthPole.end();
357  for (int idx = 0;  idx < mapper.nWeights; ++idx)
358  {
359    interpMapValue[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
360    if (iteNorthPole != interpMapValueNorthPole.find(mapper.targetWeightId[idx]))
361    {
362      interpMapValueNorthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
363    }
364
365    if (iteSouthPole != interpMapValueSouthPole.find(mapper.targetWeightId[idx]))
366    {
367      interpMapValueSouthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx]));
368    }
369  }
370  int niGloDst = domainDest_->ni_glo.getValue();
371  processPole(interpMapValueNorthPole, niGloDst);
372  processPole(interpMapValueSouthPole, niGloDst);
373
374  if (!interpMapValueNorthPole.empty())
375  {
376     std::map<int,std::vector<std::pair<int,double> > >::iterator itNorthPole = interpMapValueNorthPole.begin();
377     for (; itNorthPole != iteNorthPole; ++itNorthPole)
378     {
379       if (!(itNorthPole->second.empty()))
380        itNorthPole->second.swap(interpMapValue[itNorthPole->first]);
381     }
382  }
383
384  if (!interpMapValueSouthPole.empty())
385  {
386     std::map<int,std::vector<std::pair<int,double> > >::iterator itSouthPole = interpMapValueSouthPole.begin();
387     for (; itSouthPole != iteSouthPole; ++itSouthPole)
388     {
389       if (!(itSouthPole->second.empty()))
390        itSouthPole->second.swap(interpMapValue[itSouthPole->first]);
391     }
392  }
393
394  if (writeToFile_ && !readFromFile_) writeRemapInfo(interpMapValue);
395//  exchangeRemapInfo(interpMapValue);
396  convertRemapInfo(interpMapValue) ;
397
398  delete [] globalSrc;
399  delete [] globalSrcUnmasked;
400  delete [] globalDst;
401  delete [] globalDstUnmasked;
402
403}
404
405void CDomainAlgorithmInterpolate::processPole(std::map<int,std::vector<std::pair<int,double> > >& interMapValuePole,
406                                              int nbGlobalPointOnPole)
407{
408  CContext* context = CContext::getCurrent();
409  CContextClient* client=context->client;
410
411  ep_lib::MPI_Comm poleComme;
412  ep_lib::MPI_Comm_split(client->intraComm, interMapValuePole.empty() ? 0 : 1, 0, &poleComme);
413  if (!poleComme.is_null())
414  {
415    int nbClientPole;
416    ep_lib::MPI_Comm_size(poleComme, &nbClientPole);
417
418    std::map<int,std::vector<std::pair<int,double> > >::iterator itePole = interMapValuePole.end(), itPole,
419                                                                 itbPole = interMapValuePole.begin();
420
421    int nbWeight = 0;
422    for (itPole = itbPole; itPole != itePole; ++itPole)
423       nbWeight += itPole->second.size();
424
425    std::vector<int> recvCount(nbClientPole,0);
426    std::vector<int> displ(nbClientPole,0);
427    ep_lib::MPI_Allgather(&nbWeight,1,MPI_INT,&recvCount[0],1,MPI_INT,poleComme) ;
428
429    displ[0]=0;
430    for(int n=1;n<nbClientPole;++n) displ[n]=displ[n-1]+recvCount[n-1] ;
431    int recvSize=displ[nbClientPole-1]+recvCount[nbClientPole-1] ;
432
433    std::vector<int> sendSourceIndexBuff(nbWeight);
434    std::vector<double> sendSourceWeightBuff(nbWeight);
435    int k = 0;
436    for (itPole = itbPole; itPole != itePole; ++itPole)
437    {
438      for (int idx = 0; idx < itPole->second.size(); ++idx)
439      {
440        sendSourceIndexBuff[k] = (itPole->second)[idx].first;
441        sendSourceWeightBuff[k] = (itPole->second)[idx].second;
442        ++k;
443      }
444    }
445
446    std::vector<int> recvSourceIndexBuff(recvSize);
447    std::vector<double> recvSourceWeightBuff(recvSize);
448
449    // Gather all index and weight for pole
450    ep_lib::MPI_Allgatherv(&sendSourceIndexBuff[0],nbWeight,MPI_INT,&recvSourceIndexBuff[0],&recvCount[0],&displ[0],MPI_INT,poleComme);
451    ep_lib::MPI_Allgatherv(&sendSourceWeightBuff[0],nbWeight,MPI_DOUBLE,&recvSourceWeightBuff[0],&recvCount[0],&displ[0],MPI_DOUBLE,poleComme);
452
453    std::map<int,double> recvTemp;
454    for (int idx = 0; idx < recvSize; ++idx)
455    {
456      if (recvTemp.end() != recvTemp.find(recvSourceIndexBuff[idx]))
457        recvTemp[recvSourceIndexBuff[idx]] += recvSourceWeightBuff[idx]/nbGlobalPointOnPole;
458      else
459        recvTemp[recvSourceIndexBuff[idx]] = 0.0;
460    }
461
462    std::map<int,double>::const_iterator itRecvTemp, itbRecvTemp = recvTemp.begin(), iteRecvTemp = recvTemp.end();
463
464    for (itPole = itbPole; itPole != itePole; ++itPole)
465    {
466      itPole->second.clear();
467      for (itRecvTemp = itbRecvTemp; itRecvTemp != iteRecvTemp; ++itRecvTemp)
468          itPole->second.push_back(make_pair(itRecvTemp->first, itRecvTemp->second));
469    }
470  }
471
472}
473
474/*!
475  Compute the index mapping between domain on grid source and one on grid destination
476*/
477void CDomainAlgorithmInterpolate::computeIndexSourceMapping_(const std::vector<CArray<double,1>* >& dataAuxInputs)
478{
479  if (readFromFile_ && !writeToFile_) 
480    readRemapInfo();
481  else
482  {
483    computeRemap(); 
484  }
485}
486
487void CDomainAlgorithmInterpolate::writeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
488{ 
489  writeInterpolationInfo(fileToReadWrite_, interpMapValue);
490}
491
492void CDomainAlgorithmInterpolate::readRemapInfo()
493{ 
494  std::map<int,std::vector<std::pair<int,double> > > interpMapValue;
495  readInterpolationInfo(fileToReadWrite_, interpMapValue);
496
497  exchangeRemapInfo(interpMapValue);
498}
499
500void CDomainAlgorithmInterpolate::convertRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
501{
502  CContext* context = CContext::getCurrent();
503  CContextClient* client=context->client;
504  int clientRank = client->clientRank;
505
506  this->transformationMapping_.resize(1);
507  this->transformationWeight_.resize(1);
508
509  TransformationIndexMap& transMap = this->transformationMapping_[0];
510  TransformationWeightMap& transWeight = this->transformationWeight_[0];
511
512  std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it,
513                                                                     ite = interpMapValue.end();
514 
515  for (it = itb; it != ite; ++it)
516  {   
517    const std::vector<std::pair<int,double> >& tmp = it->second;
518    for (int i = 0; i < tmp.size(); ++i)
519    {
520      transMap[it->first].push_back(tmp[i].first);
521      transWeight[it->first].push_back(tmp[i].second);
522    }     
523  }     
524}
525
526/*!
527  Read remap information from file then distribute it among clients
528*/
529void CDomainAlgorithmInterpolate::exchangeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
530{
531  CContext* context = CContext::getCurrent();
532  CContextClient* client=context->client;
533  int clientRank = client->clientRank;
534
535  this->transformationMapping_.resize(1);
536  this->transformationWeight_.resize(1);
537
538  TransformationIndexMap& transMap = this->transformationMapping_[0];
539  TransformationWeightMap& transWeight = this->transformationWeight_[0];
540
541  boost::unordered_map<size_t,int> globalIndexOfDomainDest;
542  int ni = domainDest_->ni.getValue();
543  int nj = domainDest_->nj.getValue();
544  int ni_glo = domainDest_->ni_glo.getValue();
545  size_t globalIndex;
546  int nIndexSize = domainDest_->i_index.numElements(), i_ind, j_ind;
547  for (int idx = 0; idx < nIndexSize; ++idx)
548  {
549    i_ind=domainDest_->i_index(idx) ;
550    j_ind=domainDest_->j_index(idx) ;
551
552    globalIndex = i_ind + j_ind * ni_glo;
553    globalIndexOfDomainDest[globalIndex] = clientRank;
554  }
555
556  CClientServerMappingDistributed domainIndexClientClientMapping(globalIndexOfDomainDest,
557                                                                 client->intraComm,
558                                                                 true);
559  CArray<size_t,1> globalIndexInterp(interpMapValue.size());
560  std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it,
561                                                                     ite = interpMapValue.end();
562  size_t globalIndexCount = 0;
563  for (it = itb; it != ite; ++it)
564  {
565    globalIndexInterp(globalIndexCount) = it->first;
566    ++globalIndexCount;
567  }
568
569  domainIndexClientClientMapping.computeServerIndexMapping(globalIndexInterp);
570  const CClientServerMapping::GlobalIndexMap& globalIndexInterpSendToClient = domainIndexClientClientMapping.getGlobalIndexOnServer();
571
572  //Inform each client number of index they will receive
573  int nbClient = client->clientSize;
574  int* sendBuff = new int[nbClient];
575  int* recvBuff = new int[nbClient];
576  for (int i = 0; i < nbClient; ++i)
577  {
578    sendBuff[i] = 0;
579    recvBuff[i] = 0;
580  }
581  int sendBuffSize = 0;
582  CClientServerMapping::GlobalIndexMap::const_iterator itbMap = globalIndexInterpSendToClient.begin(), itMap,
583                                                       iteMap = globalIndexInterpSendToClient.end();
584  for (itMap = itbMap; itMap != iteMap; ++itMap)
585  {
586    const std::vector<size_t>& tmp = itMap->second;
587    int sizeIndex = 0, mapSize = (itMap->second).size();
588    for (int idx = 0; idx < mapSize; ++idx)
589    {
590//      sizeIndex += interpMapValue.at((itMap->second)[idx]).size();
591      sizeIndex += (interpMapValue[(int)(itMap->second)[idx]]).size();
592    }
593    sendBuff[itMap->first] = sizeIndex;
594    sendBuffSize += sizeIndex;
595  }
596
597
598  ep_lib::MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm);
599
600  int* sendIndexDestBuff = new int [sendBuffSize];
601  int* sendIndexSrcBuff  = new int [sendBuffSize];
602  double* sendWeightBuff = new double [sendBuffSize];
603
604  std::vector<ep_lib::MPI_Request> sendRequest(3*globalIndexInterpSendToClient.size());
605
606  int sendOffSet = 0, l = 0;
607  int position = 0;
608  for (itMap = itbMap; itMap != iteMap; ++itMap)
609  {
610    const std::vector<size_t>& indexToSend = itMap->second;
611    int mapSize = indexToSend.size();
612    int k = 0;
613    for (int idx = 0; idx < mapSize; ++idx)
614    {
615      std::vector<std::pair<int,double> >& interpMap = interpMapValue[(int)indexToSend[idx]]; //interpMapValue.at(indexToSend[idx]);
616      for (int i = 0; i < interpMap.size(); ++i)
617      {
618        sendIndexDestBuff[l] = indexToSend[idx];
619        sendIndexSrcBuff[l]  = interpMap[i].first;
620        sendWeightBuff[l]    = interpMap[i].second;
621        ++k;
622        ++l;
623      }
624    }
625
626    ep_lib::MPI_Isend(sendIndexDestBuff + sendOffSet,
627             k,
628             MPI_INT,
629             itMap->first,
630             MPI_DOMAIN_INTERPOLATION_DEST_INDEX,
631             client->intraComm,
632             &sendRequest[position++]);
633    ep_lib::MPI_Isend(sendIndexSrcBuff + sendOffSet,
634             k,
635             MPI_INT,
636             itMap->first,
637             MPI_DOMAIN_INTERPOLATION_SRC_INDEX,
638             client->intraComm,
639             &sendRequest[position++]);
640    ep_lib::MPI_Isend(sendWeightBuff + sendOffSet,
641             k,
642             MPI_DOUBLE,
643             itMap->first,
644             MPI_DOMAIN_INTERPOLATION_WEIGHT,
645             client->intraComm,
646             &sendRequest[position++]);
647    sendOffSet += k;
648  }
649
650  int recvBuffSize = recvBuff[clientRank];
651  int* recvIndexDestBuff = new int [recvBuffSize];
652  int* recvIndexSrcBuff  = new int [recvBuffSize];
653  double* recvWeightBuff = new double [recvBuffSize];
654  int receivedSize = 0;
655  int clientSrcRank;
656  while (receivedSize < recvBuffSize)
657  {
658    ep_lib::MPI_Status recvStatus;
659    ep_lib::MPI_Recv((recvIndexDestBuff + receivedSize),
660             recvBuffSize,
661             MPI_INT,
662             -2,
663             MPI_DOMAIN_INTERPOLATION_DEST_INDEX,
664             client->intraComm,
665             &recvStatus);
666
667    int countBuff = 0;
668    ep_lib::MPI_Get_count(&recvStatus, MPI_INT, &countBuff);
669    #ifdef _usingMPI
670    clientSrcRank = recvStatus.MPI_SOURCE;
671    #elif _usingEP
672    clientSrcRank = recvStatus.ep_src;
673    #endif
674
675    ep_lib::MPI_Recv((recvIndexSrcBuff + receivedSize),
676             recvBuffSize,
677             MPI_INT,
678             clientSrcRank,
679             MPI_DOMAIN_INTERPOLATION_SRC_INDEX,
680             client->intraComm,
681             &recvStatus);
682
683    ep_lib::MPI_Recv((recvWeightBuff + receivedSize),
684             recvBuffSize,
685             MPI_DOUBLE,
686             clientSrcRank,
687             MPI_DOMAIN_INTERPOLATION_WEIGHT,
688             client->intraComm,
689             &recvStatus);
690
691    for (int idx = 0; idx < countBuff; ++idx)
692    {
693      transMap[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvIndexSrcBuff + receivedSize + idx));
694      transWeight[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvWeightBuff + receivedSize + idx));
695    }
696    receivedSize += countBuff;
697  }
698
699  std::vector<ep_lib::MPI_Status> requestStatus(sendRequest.size());
700  ep_lib::MPI_Waitall(sendRequest.size(), &sendRequest[0], &requestStatus[0]);
701
702  delete [] sendIndexDestBuff;
703  delete [] sendIndexSrcBuff;
704  delete [] sendWeightBuff;
705  delete [] recvIndexDestBuff;
706  delete [] recvIndexSrcBuff;
707  delete [] recvWeightBuff;
708  delete [] sendBuff;
709  delete [] recvBuff;
710}
711 
712/*! Redefined some functions of CONetCDF4 to make use of them */
713CDomainAlgorithmInterpolate::WriteNetCdf::WriteNetCdf(const StdString& filename, const ep_lib::MPI_Comm comm)
714  : CNc4DataOutput(NULL, filename, false, false, true, comm, false, true) {}
715int CDomainAlgorithmInterpolate::WriteNetCdf::addDimensionWrite(const StdString& name, 
716                                                                const StdSize size)
717{
718  return CONetCDF4::addDimension(name, size); 
719}
720
721int CDomainAlgorithmInterpolate::WriteNetCdf::addVariableWrite(const StdString& name, nc_type type,
722                                                               const std::vector<StdString>& dim)
723{
724  return CONetCDF4::addVariable(name, type, dim);
725}
726
727void CDomainAlgorithmInterpolate::WriteNetCdf::endDefinition()
728{
729  CONetCDF4::definition_end();
730}
731
732void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<int,1>& data, const StdString& name,
733                                                              bool collective, StdSize record,
734                                                              const std::vector<StdSize>* start,
735                                                              const std::vector<StdSize>* count)
736{
737  CONetCDF4::writeData<int,1>(data, name, collective, record, start, count);
738}
739
740void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<double,1>& data, const StdString& name,
741                                                              bool collective, StdSize record,
742                                                              const std::vector<StdSize>* start,
743                                                              const std::vector<StdSize>* count)
744{
745  CONetCDF4::writeData<double,1>(data, name, collective, record, start, count);
746}
747
748/*
749   Write interpolation weights into a file
750   \param [in] filename name of output file
751   \param interpMapValue mapping of global index of domain destination and domain source as well as the corresponding weight
752*/
753void CDomainAlgorithmInterpolate::writeInterpolationInfo(std::string& filename,
754                                                         std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
755{
756  CContext* context = CContext::getCurrent();
757  CContextClient* client=context->client;
758
759  size_t n_src = domainSrc_->ni_glo * domainSrc_->nj_glo;
760  size_t n_dst = domainDest_->ni_glo * domainDest_->nj_glo;
761
762  long localNbWeight = 0;
763  long globalNbWeight;
764  long startIndex;
765  typedef std::map<int,std::vector<std::pair<int,double> > > IndexRemap;
766  IndexRemap::iterator itb = interpMapValue.begin(), it,
767                       ite = interpMapValue.end();
768  for (it = itb; it!=ite; ++it)
769  {
770    localNbWeight += (it->second).size();
771  }
772
773  CArray<int,1> src_idx(localNbWeight);
774  CArray<int,1> dst_idx(localNbWeight);
775  CArray<double,1> weights(localNbWeight);
776
777  int index = 0;
778  for (it = itb; it !=ite; ++it)
779  {
780    std::vector<std::pair<int,double> >& tmp = it->second;
781    for (int idx = 0; idx < tmp.size(); ++idx)
782    {
783      dst_idx(index) = it->first + 1;
784      src_idx(index) = tmp[idx].first + 1;
785      weights(index) = tmp[idx].second;
786      ++index;
787    }   
788  }
789
790  ep_lib::MPI_Allreduce(&localNbWeight, &globalNbWeight, 1, MPI_LONG, MPI_SUM, client->intraComm);
791  ep_lib::MPI_Scan(&localNbWeight, &startIndex, 1, MPI_LONG, MPI_SUM, client->intraComm);
792 
793  if (0 == globalNbWeight)
794  {
795    info << "There is no interpolation weights calculated between "
796         << "domain source: " << domainSrc_->getDomainOutputName()
797         << " and domain destination: " << domainDest_->getDomainOutputName()
798         << std::endl;
799    return;
800  }
801
802  std::vector<StdSize> start(1, startIndex - localNbWeight);
803  std::vector<StdSize> count(1, localNbWeight);
804 
805  WriteNetCdf netCdfWriter(filename, client->intraComm); 
806
807  // Define some dimensions
808  netCdfWriter.addDimensionWrite("n_src", n_src);
809  netCdfWriter.addDimensionWrite("n_dst", n_dst);
810  netCdfWriter.addDimensionWrite("n_weight", globalNbWeight);
811 
812  std::vector<StdString> dims(1,"n_weight");
813
814  // Add some variables
815  netCdfWriter.addVariableWrite("src_idx", NC_INT, dims);
816  netCdfWriter.addVariableWrite("dst_idx", NC_INT, dims);
817  netCdfWriter.addVariableWrite("weight", NC_DOUBLE, dims);
818
819  // End of definition
820  netCdfWriter.endDefinition();
821
822  // // Write variables
823  if (0 != localNbWeight)
824  {
825    netCdfWriter.writeDataIndex(src_idx, "src_idx", false, 0, &start, &count);
826    netCdfWriter.writeDataIndex(dst_idx, "dst_idx", false, 0, &start, &count);
827    netCdfWriter.writeDataIndex(weights, "weight", false, 0, &start, &count);
828  }
829
830  netCdfWriter.closeFile();
831}
832
833/*!
834  Read interpolation information from a file
835  \param [in] filename interpolation file
836  \param [in/out] interpMapValue Mapping between (global) index of domain on grid destination and
837         corresponding global index of domain and associated weight value on grid source
838*/
839void CDomainAlgorithmInterpolate::readInterpolationInfo(std::string& filename,
840                                                        std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
841{
842  int ncid ;
843  int weightDimId ;
844  size_t nbWeightGlo ;
845
846  CContext* context = CContext::getCurrent();
847  CContextClient* client=context->client;
848  int clientRank = client->clientRank;
849  int clientSize = client->clientSize;
850
851  nc_open(filename.c_str(),NC_NOWRITE, &ncid) ;
852  nc_inq_dimid(ncid,"n_weight",&weightDimId) ;
853  nc_inq_dimlen(ncid,weightDimId,&nbWeightGlo) ;
854
855  size_t nbWeight ;
856  size_t start ;
857  size_t div = nbWeightGlo/clientSize ;
858  size_t mod = nbWeightGlo%clientSize ;
859  if (clientRank < mod)
860  {
861    nbWeight=div+1 ;
862    start=clientRank*(div+1) ;
863  }
864  else
865  {
866    nbWeight=div ;
867    start= mod * (div+1) + (clientRank-mod) * div ;
868  }
869
870  double* weight=new double[nbWeight] ;
871  int weightId ;
872  nc_inq_varid (ncid, "weight", &weightId) ;
873  nc_get_vara_double(ncid, weightId, &start, &nbWeight, weight) ;
874
875  long* srcIndex=new long[nbWeight] ;
876  int srcIndexId ;
877  nc_inq_varid (ncid, "src_idx", &srcIndexId) ;
878  nc_get_vara_long(ncid, srcIndexId, &start, &nbWeight, srcIndex) ;
879
880  long* dstIndex=new long[nbWeight] ;
881  int dstIndexId ;
882  nc_inq_varid (ncid, "dst_idx", &dstIndexId) ;
883  nc_get_vara_long(ncid, dstIndexId, &start, &nbWeight, dstIndex) ;
884
885  for(size_t ind=0; ind<nbWeight;++ind)
886    interpMapValue[dstIndex[ind]-1].push_back(make_pair(srcIndex[ind]-1,weight[ind]));
887}
888
889}
Note: See TracBrowser for help on using the repository browser.