source: XIOS/dev/dev_olga/src/transformation/domain_algorithm_interpolate.cpp @ 1620

Last change on this file since 1620 was 1612, checked in by oabramkina, 5 years ago

Dev: adding exception handling.

To activate it, compilation flag -DXIOS_EXCEPTION should be added.

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