source: XIOS/trunk/src/transformation/domain_algorithm_interpolate.cpp @ 1615

Last change on this file since 1615 was 1615, checked in by ymipsl, 5 years ago

Interpolation enhancement :
Domain area can be used to improved global conservation when there is a discrepency between model area and xios computed area.

New domain attribute :

radius : radius af the spherical domain, used to compute area ond the sphere with a normalized radius of 1 (for remapper).

New domain_interpolate attribute :
use_area : remapper will take model computed area in order to perform a global conservation for flux integrated over the cell (mass for example).
In this case the domain attributes "area" and "radius" must be defined on the source or target domain to be taking into account.

YM

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