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

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

dev_omp

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