source: XIOS/dev/dev_trunk_omp/src/transformation/domain_algorithm_interpolate.cpp @ 1646

Last change on this file since 1646 was 1646, checked in by yushan, 5 years ago

branch merged with trunk @1645. arch file (ep&mpi) added for ADA

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