source: XIOS/dev/branch_yushan_merged/src/transformation/domain_algorithm_interpolate.cpp @ 1155

Last change on this file since 1155 was 1155, checked in by yushan, 7 years ago

test_remap OK with openmp

File size: 31.1 KB
Line 
1/*!
2   \file domain_algorithm_interpolate_from_file.cpp
3   \author Ha NGUYEN
4   \since 09 Jul 2015
5   \date 15 Sep 2015
6
7   \brief Algorithm for interpolation on a domain.
8 */
9#include "domain_algorithm_interpolate.hpp"
10#include <boost/unordered_map.hpp>
11#include "context.hpp"
12#include "context_client.hpp"
13#include "distribution_client.hpp"
14#include "client_server_mapping_distributed.hpp"
15#include "netcdf.hpp"
16#include "mapper.hpp"
17#include "mpi_tag.hpp"
18#include "domain.hpp"
19#include "grid_transformation_factory_impl.hpp"
20#include "interpolate_domain.hpp"
21#include "grid.hpp"
22
23namespace xios {
24CGenericAlgorithmTransformation* CDomainAlgorithmInterpolate::create(CGrid* gridDst, CGrid* gridSrc,
25                                                                     CTransformation<CDomain>* transformation,
26                                                                     int elementPositionInGrid,
27                                                                     std::map<int, int>& elementPositionInGridSrc2ScalarPosition,
28                                                                     std::map<int, int>& elementPositionInGridSrc2AxisPosition,
29                                                                     std::map<int, int>& elementPositionInGridSrc2DomainPosition,
30                                                                     std::map<int, int>& elementPositionInGridDst2ScalarPosition,
31                                                                     std::map<int, int>& elementPositionInGridDst2AxisPosition,
32                                                                     std::map<int, int>& elementPositionInGridDst2DomainPosition)
33{
34  std::vector<CDomain*> domainListDestP = gridDst->getDomains();
35  std::vector<CDomain*> domainListSrcP  = gridSrc->getDomains();
36
37  CInterpolateDomain* interpolateDomain = dynamic_cast<CInterpolateDomain*> (transformation);
38  int domainDstIndex = elementPositionInGridDst2DomainPosition[elementPositionInGrid];
39  int domainSrcIndex = elementPositionInGridSrc2DomainPosition[elementPositionInGrid];
40
41  return (new CDomainAlgorithmInterpolate(domainListDestP[domainDstIndex], domainListSrcP[domainSrcIndex], interpolateDomain));
42}
43
44bool CDomainAlgorithmInterpolate::registerTrans()
45{
46  CGridTransformationFactory<CDomain>::registerTransformation(TRANS_INTERPOLATE_DOMAIN, create);
47}
48
49CDomainAlgorithmInterpolate::CDomainAlgorithmInterpolate(CDomain* domainDestination, CDomain* domainSource, CInterpolateDomain* interpDomain)
50: CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain), writeToFile_(false), readFromFile_(false)
51{
52  CContext* context = CContext::getCurrent();
53  interpDomain_->checkValid(domainSource);
54  fileToReadWrite_ = "xios_interpolation_weights_";
55
56  if (interpDomain_->weight_filename.isEmpty())
57  {
58    fileToReadWrite_ += context->getId() + "_" + 
59                    domainSource->getDomainOutputName() + "_" + 
60                    domainDestination->getDomainOutputName() + ".nc";   
61  }
62  else 
63    fileToReadWrite_ = interpDomain_->weight_filename;
64
65  ifstream f(fileToReadWrite_.c_str()); 
66  switch (interpDomain_->mode)
67  {
68    case CInterpolateDomain::mode_attr::read:
69      readFromFile_ = true;     
70      break;
71    case CInterpolateDomain::mode_attr::compute:
72      readFromFile_ = false;
73      break;
74    case CInterpolateDomain::mode_attr::read_or_compute:     
75      if (!f.good())
76        readFromFile_ = false;
77      else
78        readFromFile_ = true;
79      break;
80    default:
81      break;
82  } 
83
84  writeToFile_ = interpDomain_->write_weight; 
85   
86}
87
88/*!
89  Compute remap with integrated remap calculation module
90*/
91void CDomainAlgorithmInterpolate::computeRemap()
92{
93  using namespace sphereRemap;
94
95  CContext* context = CContext::getCurrent();
96  CContextClient* client=context->client;
97  int clientRank = client->clientRank;
98  int i, j, k, idx;
99  std::vector<double> srcPole(3,0), dstPole(3,0);
100  int orderInterp = interpDomain_->order.getValue();
101  bool renormalize ;
102  bool quantity ;
103
104  if (interpDomain_->renormalize.isEmpty()) renormalize=true;
105  else renormalize = interpDomain_->renormalize;
106
107  if (interpDomain_->quantity.isEmpty()) quantity=false;
108  else quantity = interpDomain_->quantity;
109
110  const double poleValue = 90.0;
111  const int constNVertex = 4; // Value by default number of vertex for rectangular domain
112  int nVertexSrc, nVertexDest;
113  nVertexSrc = constNVertex;
114  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_)
392     writeRemapInfo(interpMapValue);
393  exchangeRemapInfo(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 split_key;
408  ep_lib::MPI_Comm_rank(client->intraComm, &split_key);
409
410  ep_lib::MPI_Comm poleComme(MPI_COMM_NULL);
411  ep_lib::MPI_Comm_split(client->intraComm, interMapValuePole.empty() ? 9 : 1, split_key, &poleComme);
412  if (MPI_COMM_NULL != poleComme)
413  {
414    int nbClientPole;
415    ep_lib::MPI_Comm_size(poleComme, &nbClientPole);
416
417    std::map<int,std::vector<std::pair<int,double> > >::iterator itePole = interMapValuePole.end(), itPole,
418                                                                 itbPole = interMapValuePole.begin();
419
420    int nbWeight = 0;
421    for (itPole = itbPole; itPole != itePole; ++itPole)
422       nbWeight += itPole->second.size();
423
424    std::vector<int> recvCount(nbClientPole,0);
425    std::vector<int> displ(nbClientPole,0);
426    ep_lib::MPI_Allgather(&nbWeight,1,MPI_INT,&recvCount[0],1,MPI_INT,poleComme) ;
427
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
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
500
501/*!
502  Read remap information from file then distribute it among clients
503*/
504void CDomainAlgorithmInterpolate::exchangeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
505{
506  CContext* context = CContext::getCurrent();
507  CContextClient* client=context->client;
508  int clientRank = client->clientRank;
509
510  this->transformationMapping_.resize(1);
511  this->transformationWeight_.resize(1);
512
513  TransformationIndexMap& transMap = this->transformationMapping_[0];
514  TransformationWeightMap& transWeight = this->transformationWeight_[0];
515
516  boost::unordered_map<size_t,int> globalIndexOfDomainDest;
517  int ni = domainDest_->ni.getValue();
518  int nj = domainDest_->nj.getValue();
519  int ni_glo = domainDest_->ni_glo.getValue();
520  size_t globalIndex;
521  int nIndexSize = domainDest_->i_index.numElements(), i_ind, j_ind;
522  for (int idx = 0; idx < nIndexSize; ++idx)
523  {
524    i_ind=domainDest_->i_index(idx) ;
525    j_ind=domainDest_->j_index(idx) ;
526
527    globalIndex = i_ind + j_ind * ni_glo;
528    globalIndexOfDomainDest[globalIndex] = clientRank;
529  }
530
531  CClientServerMappingDistributed domainIndexClientClientMapping(globalIndexOfDomainDest,
532                                                                 client->intraComm,
533                                                                 true);
534  CArray<size_t,1> globalIndexInterp(interpMapValue.size());
535  std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it,
536                                                                     ite = interpMapValue.end();
537  size_t globalIndexCount = 0;
538  for (it = itb; it != ite; ++it)
539  {
540    globalIndexInterp(globalIndexCount) = it->first;
541    ++globalIndexCount;
542  }
543
544  domainIndexClientClientMapping.computeServerIndexMapping(globalIndexInterp);
545  const CClientServerMapping::GlobalIndexMap& globalIndexInterpSendToClient = domainIndexClientClientMapping.getGlobalIndexOnServer();
546
547  //Inform each client number of index they will receive
548  int nbClient = client->clientSize;
549  int* sendBuff = new int[nbClient];
550  int* recvBuff = new int[nbClient];
551  for (int i = 0; i < nbClient; ++i)
552  {
553    sendBuff[i] = 0;
554    recvBuff[i] = 0;
555  }
556  int sendBuffSize = 0;
557  CClientServerMapping::GlobalIndexMap::const_iterator itbMap = globalIndexInterpSendToClient.begin(), itMap,
558                                                       iteMap = globalIndexInterpSendToClient.end();
559  for (itMap = itbMap; itMap != iteMap; ++itMap)
560  {
561    const std::vector<size_t>& tmp = itMap->second;
562    int sizeIndex = 0, mapSize = (itMap->second).size();
563    for (int idx = 0; idx < mapSize; ++idx)
564    {
565//      sizeIndex += interpMapValue.at((itMap->second)[idx]).size();
566      sizeIndex += (interpMapValue[(int)(itMap->second)[idx]]).size();
567    }
568    sendBuff[itMap->first] = sizeIndex;
569    sendBuffSize += sizeIndex;
570  }
571
572
573  ep_lib::MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm);
574
575  int* sendIndexDestBuff = new int [sendBuffSize];
576  int* sendIndexSrcBuff  = new int [sendBuffSize];
577  double* sendWeightBuff = new double [sendBuffSize];
578
579  std::vector<ep_lib::MPI_Request> sendRequest;
580
581  int sendOffSet = 0, l = 0;
582  for (itMap = itbMap; itMap != iteMap; ++itMap)
583  {
584    const std::vector<size_t>& indexToSend = itMap->second;
585    int mapSize = indexToSend.size();
586    int k = 0;
587    for (int idx = 0; idx < mapSize; ++idx)
588    {
589      std::vector<std::pair<int,double> >& interpMap = interpMapValue[(int)indexToSend[idx]]; //interpMapValue.at(indexToSend[idx]);
590      for (int i = 0; i < interpMap.size(); ++i)
591      {
592        sendIndexDestBuff[l] = indexToSend[idx];
593        sendIndexSrcBuff[l]  = interpMap[i].first;
594        sendWeightBuff[l]    = interpMap[i].second;
595        ++k;
596        ++l;
597      }
598    }
599
600    sendRequest.push_back(ep_lib::MPI_Request());
601    ep_lib::MPI_Isend(sendIndexDestBuff + sendOffSet,
602             k,
603             MPI_INT,
604             itMap->first,
605             MPI_DOMAIN_INTERPOLATION_DEST_INDEX,
606             client->intraComm,
607             &sendRequest.back());
608    sendRequest.push_back(ep_lib::MPI_Request());
609    ep_lib::MPI_Isend(sendIndexSrcBuff + sendOffSet,
610             k,
611             MPI_INT,
612             itMap->first,
613             MPI_DOMAIN_INTERPOLATION_SRC_INDEX,
614             client->intraComm,
615             &sendRequest.back());
616    sendRequest.push_back(ep_lib::MPI_Request());
617    ep_lib::MPI_Isend(sendWeightBuff + sendOffSet,
618             k,
619             MPI_DOUBLE,
620             itMap->first,
621             MPI_DOMAIN_INTERPOLATION_WEIGHT,
622             client->intraComm,
623             &sendRequest.back());
624    sendOffSet += k;
625  }
626
627  int recvBuffSize = recvBuff[clientRank];
628  int* recvIndexDestBuff = new int [recvBuffSize];
629  int* recvIndexSrcBuff  = new int [recvBuffSize];
630  double* recvWeightBuff = new double [recvBuffSize];
631  int receivedSize = 0;
632  int clientSrcRank;
633  while (receivedSize < recvBuffSize)
634  {
635    ep_lib::MPI_Status recvStatus;
636    ep_lib::MPI_Recv((recvIndexDestBuff + receivedSize),
637             recvBuffSize,
638             MPI_INT,
639             MPI_ANY_SOURCE,
640             MPI_DOMAIN_INTERPOLATION_DEST_INDEX,
641             client->intraComm,
642             &recvStatus);
643
644    int countBuff = 0;
645    ep_lib::MPI_Get_count(&recvStatus, MPI_INT, &countBuff);
646    #ifdef _usingMPI
647    clientSrcRank = recvStatus.MPI_SOURCE;
648    #elif _usingEP
649    clientSrcRank = recvStatus.ep_src;
650    #endif
651    ep_lib::MPI_Recv((recvIndexSrcBuff + receivedSize),
652             recvBuffSize,
653             MPI_INT,
654             clientSrcRank,
655             MPI_DOMAIN_INTERPOLATION_SRC_INDEX,
656             client->intraComm,
657             &recvStatus);
658
659    ep_lib::MPI_Recv((recvWeightBuff + receivedSize),
660             recvBuffSize,
661             MPI_DOUBLE,
662             clientSrcRank,
663             MPI_DOMAIN_INTERPOLATION_WEIGHT,
664             client->intraComm,
665             &recvStatus);
666
667    for (int idx = 0; idx < countBuff; ++idx)
668    {
669      transMap[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvIndexSrcBuff + receivedSize + idx));
670      transWeight[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvWeightBuff + receivedSize + idx));
671    }
672    receivedSize += countBuff;
673  }
674
675  std::vector<ep_lib::MPI_Status> requestStatus(sendRequest.size());
676  ep_lib::MPI_Status stat_ignore;
677  ep_lib::MPI_Waitall(sendRequest.size(), &sendRequest[0], &stat_ignore);
678
679  delete [] sendIndexDestBuff;
680  delete [] sendIndexSrcBuff;
681  delete [] sendWeightBuff;
682  delete [] recvIndexDestBuff;
683  delete [] recvIndexSrcBuff;
684  delete [] recvWeightBuff;
685  delete [] sendBuff;
686  delete [] recvBuff;
687}
688 
689/*! Redefined some functions of CONetCDF4 to make use of them */
690CDomainAlgorithmInterpolate::WriteNetCdf::WriteNetCdf(const StdString& filename, const MPI_Comm comm)
691  : CNc4DataOutput(NULL, filename, false, false, true, comm, false, true) {}
692int CDomainAlgorithmInterpolate::WriteNetCdf::addDimensionWrite(const StdString& name, 
693                                                                const StdSize size)
694{
695  return CONetCDF4::addDimension(name, size); 
696}
697
698int CDomainAlgorithmInterpolate::WriteNetCdf::addVariableWrite(const StdString& name, nc_type type,
699                                                               const std::vector<StdString>& dim)
700{
701  return CONetCDF4::addVariable(name, type, dim);
702}
703
704void CDomainAlgorithmInterpolate::WriteNetCdf::endDefinition()
705{
706  CONetCDF4::definition_end();
707}
708
709void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<int,1>& data, const StdString& name,
710                                                              bool collective, StdSize record,
711                                                              const std::vector<StdSize>* start,
712                                                              const std::vector<StdSize>* count)
713{
714  CONetCDF4::writeData<int,1>(data, name, collective, record, start, count);
715}
716
717void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<double,1>& data, const StdString& name,
718                                                              bool collective, StdSize record,
719                                                              const std::vector<StdSize>* start,
720                                                              const std::vector<StdSize>* count)
721{
722  CONetCDF4::writeData<double,1>(data, name, collective, record, start, count);
723}
724
725/*
726   Write interpolation weights into a file
727   \param [in] filename name of output file
728   \param interpMapValue mapping of global index of domain destination and domain source as well as the corresponding weight
729*/
730void CDomainAlgorithmInterpolate::writeInterpolationInfo(std::string& filename,
731                                                         std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
732{
733  CContext* context = CContext::getCurrent();
734  CContextClient* client=context->client;
735
736  size_t n_src = domainSrc_->ni_glo * domainSrc_->nj_glo;
737  size_t n_dst = domainDest_->ni_glo * domainDest_->nj_glo;
738
739  long localNbWeight = 0;
740  long globalNbWeight;
741  long startIndex;
742  typedef std::map<int,std::vector<std::pair<int,double> > > IndexRemap;
743  IndexRemap::iterator itb = interpMapValue.begin(), it,
744                       ite = interpMapValue.end();
745  for (it = itb; it!=ite; ++it)
746  {
747    localNbWeight += (it->second).size();
748  }
749
750  CArray<int,1> src_idx(localNbWeight);
751  CArray<int,1> dst_idx(localNbWeight);
752  CArray<double,1> weights(localNbWeight);
753
754  int index = 0;
755  for (it = itb; it !=ite; ++it)
756  {
757    std::vector<std::pair<int,double> >& tmp = it->second;
758    for (int idx = 0; idx < tmp.size(); ++idx)
759    {
760      dst_idx(index) = it->first + 1;
761      src_idx(index) = tmp[idx].first + 1;
762      weights(index) = tmp[idx].second;
763      ++index;
764    }   
765  }
766
767  ep_lib::MPI_Allreduce(&localNbWeight, &globalNbWeight, 1, MPI_LONG, MPI_SUM, client->intraComm);
768  ep_lib::MPI_Scan(&localNbWeight, &startIndex, 1, MPI_LONG, MPI_SUM, client->intraComm);
769 
770  if (0 == globalNbWeight)
771  {
772    info << "There is no interpolation weights calculated between "
773         << "domain source: " << domainSrc_->getDomainOutputName()
774         << " and domain destination: " << domainDest_->getDomainOutputName()
775         << std::endl;
776    return;
777  }
778
779  std::vector<StdSize> start(1, startIndex - localNbWeight);
780  std::vector<StdSize> count(1, localNbWeight);
781
782  WriteNetCdf netCdfWriter(filename, static_cast<MPI_Comm>(client->intraComm.mpi_comm));
783
784  // Define some dimensions
785  netCdfWriter.addDimensionWrite("n_src", n_src);
786  netCdfWriter.addDimensionWrite("n_dst", n_dst);
787  netCdfWriter.addDimensionWrite("n_weight", globalNbWeight);
788 
789  std::vector<StdString> dims(1,"n_weight");
790
791  // Add some variables
792  netCdfWriter.addVariableWrite("src_idx", NC_INT, dims);
793  netCdfWriter.addVariableWrite("dst_idx", NC_INT, dims);
794  netCdfWriter.addVariableWrite("weight", NC_DOUBLE, dims);
795
796  // End of definition
797  netCdfWriter.endDefinition();
798
799  // // Write variables
800  if (0 != localNbWeight)
801  {
802    netCdfWriter.writeDataIndex(src_idx, "src_idx", false, 0, &start, &count);
803    netCdfWriter.writeDataIndex(dst_idx, "dst_idx", false, 0, &start, &count);
804    netCdfWriter.writeDataIndex(weights, "weight", false, 0, &start, &count);
805  }
806
807  netCdfWriter.closeFile();
808}
809
810/*!
811  Read interpolation information from a file
812  \param [in] filename interpolation file
813  \param [in/out] interpMapValue Mapping between (global) index of domain on grid destination and
814         corresponding global index of domain and associated weight value on grid source
815*/
816void CDomainAlgorithmInterpolate::readInterpolationInfo(std::string& filename,
817                                                        std::map<int,std::vector<std::pair<int,double> > >& interpMapValue)
818{
819  int ncid ;
820  int weightDimId ;
821  size_t nbWeightGlo ;
822
823  CContext* context = CContext::getCurrent();
824  CContextClient* client=context->client;
825  int clientRank = client->clientRank;
826  int clientSize = client->clientSize;
827
828  nc_open(filename.c_str(),NC_NOWRITE, &ncid) ;
829  nc_inq_dimid(ncid,"n_weight",&weightDimId) ;
830  nc_inq_dimlen(ncid,weightDimId,&nbWeightGlo) ;
831
832  size_t nbWeight ;
833  size_t start ;
834  size_t div = nbWeightGlo/clientSize ;
835  size_t mod = nbWeightGlo%clientSize ;
836  if (clientRank < mod)
837  {
838    nbWeight=div+1 ;
839    start=clientRank*(div+1) ;
840  }
841  else
842  {
843    nbWeight=div ;
844    start= mod * (div+1) + (clientRank-mod) * div ;
845  }
846
847  double* weight=new double[nbWeight] ;
848  int weightId ;
849  nc_inq_varid (ncid, "weight", &weightId) ;
850  nc_get_vara_double(ncid, weightId, &start, &nbWeight, weight) ;
851
852  long* srcIndex=new long[nbWeight] ;
853  int srcIndexId ;
854  nc_inq_varid (ncid, "src_idx", &srcIndexId) ;
855  nc_get_vara_long(ncid, srcIndexId, &start, &nbWeight, srcIndex) ;
856
857  long* dstIndex=new long[nbWeight] ;
858  int dstIndexId ;
859  nc_inq_varid (ncid, "dst_idx", &dstIndexId) ;
860  nc_get_vara_long(ncid, dstIndexId, &start, &nbWeight, dstIndex) ;
861
862  for(size_t ind=0; ind<nbWeight;++ind)
863    interpMapValue[dstIndex[ind]-1].push_back(make_pair(srcIndex[ind]-1,weight[ind]));
864}
865
866}
Note: See TracBrowser for help on using the repository browser.