source: XIOS/trunk/src/node/domain.cpp @ 1578

Last change on this file since 1578 was 1578, checked in by ymipsl, 6 years ago

Bug fix when reading rectilinear and attemp to call generate_rectilinear_domain filter. Coordinates in file was not taking into account.

YM

  • Property copyright set to
    Software name : XIOS (Xml I/O Server)
    http://forge.ipsl.jussieu.fr/ioserver
    Creation date : January 2009
    Licence : CeCCIL version2
    see license file in root directory : Licence_CeCILL_V2-en.txt
    or http://www.cecill.info/licences/Licence_CeCILL_V2-en.html
    Holder : CEA/LSCE (Laboratoire des Sciences du CLimat et de l'Environnement)
    CNRS/IPSL (Institut Pierre Simon Laplace)
    Project Manager : Yann Meurdesoif
    yann.meurdesoif@cea.fr
  • Property svn:executable set to *
File size: 109.4 KB
Line 
1#include "domain.hpp"
2
3#include "attribute_template.hpp"
4#include "object_template.hpp"
5#include "group_template.hpp"
6
7#include "xios_spl.hpp"
8#include "event_client.hpp"
9#include "event_server.hpp"
10#include "buffer_in.hpp"
11#include "message.hpp"
12#include "type.hpp"
13#include "context.hpp"
14#include "context_client.hpp"
15#include "context_server.hpp"
16#include "array_new.hpp"
17#include "distribution_client.hpp"
18#include "server_distribution_description.hpp"
19#include "client_server_mapping_distributed.hpp"
20
21#include <algorithm>
22
23namespace xios {
24
25   /// ////////////////////// Définitions ////////////////////// ///
26
27   CDomain::CDomain(void)
28      : CObjectTemplate<CDomain>(), CDomainAttributes()
29      , isChecked(false), relFiles(), isClientChecked(false), nbSenders(), indSrv_(), connectedServerRank_()
30      , hasBounds(false), hasArea(false), isCompressible_(false), isUnstructed_(false)
31      , isClientAfterTransformationChecked(false), hasLonLat(false)
32      , isRedistributed_(false), hasPole(false)
33      , lonvalue(), latvalue(), bounds_lonvalue(), bounds_latvalue()
34      , globalLocalIndexMap_(), computedWrittenIndex_(false)
35      , clients(), hasLatInReadFile_(false), hasBoundsLatInReadFile_(false)
36      , hasLonInReadFile_(false), hasBoundsLonInReadFile_(false)
37   {
38   }
39
40   CDomain::CDomain(const StdString & id)
41      : CObjectTemplate<CDomain>(id), CDomainAttributes()
42      , isChecked(false), relFiles(), isClientChecked(false), nbSenders(), indSrv_(), connectedServerRank_() 
43      , hasBounds(false), hasArea(false), isCompressible_(false), isUnstructed_(false)
44      , isClientAfterTransformationChecked(false), hasLonLat(false)
45      , isRedistributed_(false), hasPole(false)
46      , lonvalue(), latvalue(), bounds_lonvalue(), bounds_latvalue()
47      , globalLocalIndexMap_(), computedWrittenIndex_(false)
48      , clients(), hasLatInReadFile_(false), hasBoundsLatInReadFile_(false)
49      , hasLonInReadFile_(false), hasBoundsLonInReadFile_(false)
50   {
51    }
52
53   CDomain::~CDomain(void)
54   {
55   }
56
57   ///---------------------------------------------------------------
58
59   void CDomain::assignMesh(const StdString meshName, const int nvertex)
60   {
61     mesh = CMesh::getMesh(meshName, nvertex);
62   }
63
64   CDomain* CDomain::createDomain()
65   {
66     CDomain* domain = CDomainGroup::get("domain_definition")->createChild();
67     return domain;
68   }
69
70   std::map<StdString, ETranformationType> CDomain::transformationMapList_ = std::map<StdString, ETranformationType>();
71   bool CDomain::_dummyTransformationMapList = CDomain::initializeTransformationMap(CDomain::transformationMapList_);
72
73   bool CDomain::initializeTransformationMap(std::map<StdString, ETranformationType>& m)
74   {
75     m["zoom_domain"] = TRANS_ZOOM_DOMAIN;
76     m["interpolate_domain"] = TRANS_INTERPOLATE_DOMAIN;
77     m["generate_rectilinear_domain"] = TRANS_GENERATE_RECTILINEAR_DOMAIN;
78     m["compute_connectivity_domain"] = TRANS_COMPUTE_CONNECTIVITY_DOMAIN;
79     m["expand_domain"] = TRANS_EXPAND_DOMAIN;
80     m["reorder_domain"] = TRANS_REORDER_DOMAIN;
81     m["extract_domain"] = TRANS_EXTRACT_DOMAIN;
82   }
83
84   const std::set<StdString> & CDomain::getRelFiles(void) const
85   {
86      return (this->relFiles);
87   }
88
89
90   /*!
91     Returns the number of indexes written by each server.
92     \return the number of indexes written by each server
93   */
94   int CDomain::getNumberWrittenIndexes(MPI_Comm writtenCom)
95   {
96     int writtenSize;
97     MPI_Comm_size(writtenCom, &writtenSize);
98     return numberWrittenIndexes_[writtenSize];
99   }
100
101   /*!
102     Returns the total number of indexes written by the servers.
103     \return the total number of indexes written by the servers
104   */
105   int CDomain::getTotalNumberWrittenIndexes(MPI_Comm writtenCom)
106   {
107     int writtenSize;
108     MPI_Comm_size(writtenCom, &writtenSize);
109     return totalNumberWrittenIndexes_[writtenSize];
110   }
111
112   /*!
113     Returns the offset of indexes written by each server.
114     \return the offset of indexes written by each server
115   */
116   int CDomain::getOffsetWrittenIndexes(MPI_Comm writtenCom)
117   {
118     int writtenSize;
119     MPI_Comm_size(writtenCom, &writtenSize);
120     return offsetWrittenIndexes_[writtenSize];
121   }
122
123   CArray<int, 1>& CDomain::getCompressedIndexToWriteOnServer(MPI_Comm writtenCom)
124   {
125     int writtenSize;
126     MPI_Comm_size(writtenCom, &writtenSize);
127     return compressedIndexToWriteOnServer[writtenSize];
128   }
129
130   //----------------------------------------------------------------
131
132   /*!
133    * Compute the minimum buffer size required to send the attributes to the server(s).
134    *
135    * \return A map associating the server rank with its minimum buffer size.
136    */
137   std::map<int, StdSize> CDomain::getAttributesBufferSize(CContextClient* client, bool bufferForWriting /*= false*/)
138   {
139
140     std::map<int, StdSize> attributesSizes = getMinimumBufferSizeForAttributes(client);
141
142     if (client->isServerLeader())
143     {
144       // size estimation for sendDistributionAttribut
145       size_t size = 11 * sizeof(size_t);
146
147       const std::list<int>& ranks = client->getRanksServerLeader();
148       for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
149       {
150         if (size > attributesSizes[*itRank])
151           attributesSizes[*itRank] = size;
152       }
153     }
154
155     std::unordered_map<int, vector<size_t> >::const_iterator itIndexEnd = indSrv_[client->serverSize].end();
156     // std::map<int, std::vector<int> >::const_iterator itWrittenIndexEnd = indWrittenSrv_.end();
157     for (size_t k = 0; k < connectedServerRank_[client->serverSize].size(); ++k)
158     {
159       int rank = connectedServerRank_[client->serverSize][k];
160       std::unordered_map<int, std::vector<size_t> >::const_iterator it = indSrv_[client->serverSize].find(rank);
161       size_t idxCount = (it != itIndexEnd) ? it->second.size() : 0;
162
163       // size estimation for sendIndex (and sendArea which is always smaller or equal)
164       size_t sizeIndexEvent = 2 * sizeof(size_t) + 2 * CArray<int,1>::size(idxCount);
165       // if (isCompressible_)
166       // {
167       //   std::map<int, std::vector<int> >::const_iterator itWritten = indWrittenSrv_.find(rank);
168       //   size_t writtenIdxCount = (itWritten != itWrittenIndexEnd) ? itWritten->second.size() : 0;
169       //   sizeIndexEvent += CArray<int,1>::size(writtenIdxCount);
170       // }
171
172       // size estimation for sendLonLat
173       size_t sizeLonLatEvent = CArray<double,1>::size(idxCount);
174       if (hasBounds)
175         sizeLonLatEvent += CArray<double,2>::size(nvertex * idxCount);
176
177       size_t size = CEventClient::headerSize + getId().size() + sizeof(size_t) + std::max(sizeIndexEvent, sizeLonLatEvent);
178       if (size > attributesSizes[rank])
179         attributesSizes[rank] = size;
180     }
181
182     return attributesSizes;
183   }
184
185   //----------------------------------------------------------------
186
187   bool CDomain::isEmpty(void) const
188   {
189     return ((this->i_index.isEmpty()) || (0 == this->i_index.numElements()));
190   }
191
192   //----------------------------------------------------------------
193
194   bool CDomain::IsWritten(const StdString & filename) const
195   {
196      return (this->relFiles.find(filename) != this->relFiles.end());
197   }
198
199   bool CDomain::isWrittenCompressed(const StdString& filename) const
200   {
201      return (this->relFilesCompressed.find(filename) != this->relFilesCompressed.end());
202   }
203
204   //----------------------------------------------------------------
205
206   bool CDomain::isDistributed(void) const
207   {
208      bool distributed =  !((!ni.isEmpty() && (ni == ni_glo) && !nj.isEmpty() && (nj == nj_glo)) ||
209              (!i_index.isEmpty() && i_index.numElements() == ni_glo*nj_glo));
210      distributed |= (1 == CContext::getCurrent()->client->clientSize);
211
212      return distributed;
213   }
214
215   //----------------------------------------------------------------
216
217   /*!
218    * Test whether the data defined on the domain can be outputted in a compressed way.
219    *
220    * \return true if and only if a mask was defined for this domain
221    */
222   bool CDomain::isCompressible(void) const
223   {
224      return isCompressible_;
225   }
226
227   void CDomain::addRelFile(const StdString & filename)
228   {
229      this->relFiles.insert(filename);
230   }
231
232   void CDomain::addRelFileCompressed(const StdString& filename)
233   {
234      this->relFilesCompressed.insert(filename);
235   }
236
237   StdString CDomain::GetName(void)   { return (StdString("domain")); }
238   StdString CDomain::GetDefName(void){ return (CDomain::GetName()); }
239   ENodeType CDomain::GetType(void)   { return (eDomain); }
240
241   //----------------------------------------------------------------
242
243   /*!
244      Verify if all distribution information of a domain are available
245      This checking verifies the definition of distribution attributes (ni, nj, ibegin, jbegin)
246   */
247   bool CDomain::distributionAttributesHaveValue() const
248   {
249      bool hasValues = true;
250
251      if (ni.isEmpty() && ibegin.isEmpty() && i_index.isEmpty())
252      {
253        hasValues = false;
254        return hasValues;
255      }
256
257      return hasValues;
258   }
259
260   /*!
261     Redistribute RECTILINEAR or CURVILINEAR domain with a number of local domains.
262   All attributes ni,nj,ibegin,jbegin (if defined) will be rewritten
263   The optional attributes lonvalue, latvalue will be added. Because this function only serves (for now)
264   for interpolation from unstructured domain to rectilinear one, range of latvalue is 0-360 and lonvalue is -90 - +90
265    \param [in] nbLocalDomain number of local domain on the domain destination
266   */
267   void CDomain::redistribute(int nbLocalDomain)
268   {
269     if (this->isRedistributed_) return;
270
271     this->isRedistributed_ = true;
272     CContext* context = CContext::getCurrent();
273     // For now the assumption is that secondary server pools consist of the same number of procs.
274     // CHANGE the line below if the assumption changes.
275     CContextClient* client = (0 != context->clientPrimServer.size()) ? context->clientPrimServer[0] : context->client;
276     int rankClient = client->clientRank;
277     int rankOnDomain = rankClient%nbLocalDomain;
278
279     if (ni_glo.isEmpty() || ni_glo <= 0 )
280     {
281        ERROR("CDomain::redistribute(int nbLocalDomain)",
282           << "[ Id = " << this->getId() << " ] "
283           << "The global domain is badly defined,"
284           << " check the \'ni_glo\'  value !")
285     }
286
287     if (nj_glo.isEmpty() || nj_glo <= 0 )
288     {
289        ERROR("CDomain::redistribute(int nbLocalDomain)",
290           << "[ Id = " << this->getId() << " ] "
291           << "The global domain is badly defined,"
292           << " check the \'nj_glo\'  value !")
293     }
294
295     if ((type_attr::rectilinear == type)  || (type_attr::curvilinear == type))
296     {
297        int globalDomainSize = ni_glo * nj_glo;
298        if (globalDomainSize <= nbLocalDomain)
299        {
300          for (int idx = 0; idx < nbLocalDomain; ++idx)
301          {
302            if (rankOnDomain < globalDomainSize)
303            {
304              int iIdx = rankOnDomain % ni_glo;
305              int jIdx = rankOnDomain / ni_glo;
306              ibegin.setValue(iIdx); jbegin.setValue(jIdx);
307              ni.setValue(1); nj.setValue(1);
308            }
309            else
310            {
311              ibegin.setValue(0); jbegin.setValue(0);
312              ni.setValue(0); nj.setValue(0);
313            }
314          }
315        }
316        else
317        {
318          float njGlo = nj_glo.getValue();
319          float niGlo = ni_glo.getValue();
320          int nbProcOnX, nbProcOnY, range;
321
322          // Compute (approximately) number of segment on x and y axis
323          float yOverXRatio = njGlo/niGlo;
324
325          nbProcOnX = std::ceil(std::sqrt(nbLocalDomain/yOverXRatio));
326          nbProcOnY = std::ceil(((float)nbLocalDomain)/nbProcOnX);
327
328          // Simple distribution: Sweep from top to bottom, left to right
329          // Calculate local begin on x
330          std::vector<int> ibeginVec(nbProcOnX,0), jbeginVec(nbProcOnY,0);
331          std::vector<int> niVec(nbProcOnX), njVec(nbProcOnY);
332          for (int i = 1; i < nbProcOnX; ++i)
333          {
334            range = ni_glo / nbProcOnX;
335            if (i < (ni_glo%nbProcOnX)) ++range;
336            niVec[i-1] = range;
337            ibeginVec[i] = ibeginVec[i-1] + niVec[i-1];
338          }
339          niVec[nbProcOnX-1] = ni_glo - ibeginVec[nbProcOnX-1];
340
341          // Calculate local begin on y
342          for (int j = 1; j < nbProcOnY; ++j)
343          {
344            range = nj_glo / nbProcOnY;
345            if (j < (nj_glo%nbProcOnY)) ++range;
346            njVec[j-1] = range;
347            jbeginVec[j] = jbeginVec[j-1] + njVec[j-1];
348          }
349          njVec[nbProcOnY-1] = nj_glo - jbeginVec[nbProcOnY-1];
350
351          // Now assign value to ni, ibegin, nj, jbegin
352          int iIdx = rankOnDomain % nbProcOnX;
353          int jIdx = rankOnDomain / nbProcOnX;
354
355          if (rankOnDomain != (nbLocalDomain-1))
356          {
357            ibegin.setValue(ibeginVec[iIdx]);
358            jbegin.setValue(jbeginVec[jIdx]);
359            nj.setValue(njVec[jIdx]);
360            ni.setValue(niVec[iIdx]);
361          }
362          else // just merge all the remaining rectangle into the last one
363          {
364            ibegin.setValue(ibeginVec[iIdx]);
365            jbegin.setValue(jbeginVec[jIdx]);
366            nj.setValue(njVec[jIdx]);
367            ni.setValue(ni_glo - ibeginVec[iIdx]);
368          }
369        } 
370     }
371     else  // unstructured domain
372     {
373       if (this->i_index.isEmpty())
374       {
375          int globalDomainSize = ni_glo * nj_glo;
376          if (globalDomainSize <= nbLocalDomain)
377          {
378            for (int idx = 0; idx < nbLocalDomain; ++idx)
379            {
380              if (rankOnDomain < globalDomainSize)
381              {
382                int iIdx = rankOnDomain % ni_glo;
383                int jIdx = rankOnDomain / ni_glo;
384                ibegin.setValue(iIdx); jbegin.setValue(jIdx);
385                ni.setValue(1); nj.setValue(1);
386              }
387              else
388              {
389                ibegin.setValue(0); jbegin.setValue(0);
390                ni.setValue(0); nj.setValue(0);
391              }
392            }
393          }
394          else
395          {
396            float njGlo = nj_glo.getValue();
397            float niGlo = ni_glo.getValue();
398            std::vector<int> ibeginVec(nbLocalDomain,0);
399            std::vector<int> niVec(nbLocalDomain);
400            for (int i = 1; i < nbLocalDomain; ++i)
401            {
402              int range = ni_glo / nbLocalDomain;
403              if (i < (ni_glo%nbLocalDomain)) ++range;
404              niVec[i-1] = range;
405              ibeginVec[i] = ibeginVec[i-1] + niVec[i-1];
406            }
407            niVec[nbLocalDomain-1] = ni_glo - ibeginVec[nbLocalDomain-1];
408
409            int iIdx = rankOnDomain % nbLocalDomain;
410            ibegin.setValue(ibeginVec[iIdx]);
411            jbegin.setValue(0);
412            ni.setValue(niVec[iIdx]);
413            nj.setValue(1);
414          }
415
416          i_index.resize(ni);         
417          for(int idx = 0; idx < ni; ++idx) i_index(idx)=ibegin+idx;
418        }
419        else
420        {
421          ibegin.setValue(this->i_index(0));
422          jbegin.setValue(0);
423          ni.setValue(this->i_index.numElements());
424          nj.setValue(1);
425        }
426     }
427
428     checkDomain();
429   }
430
431   /*!
432     Fill in longitude and latitude whose values are read from file
433   */
434   void CDomain::fillInLonLat()
435   {
436     switch (type)
437     {
438      case type_attr::rectilinear:
439        fillInRectilinearLonLat();
440        break;
441      case type_attr::curvilinear:
442        fillInCurvilinearLonLat();
443        break;
444      case type_attr::unstructured:
445        fillInUnstructuredLonLat();
446        break;
447
448      default:
449      break;
450     }
451     completeLonLatClient() ;
452
453   }
454
455   /*!
456     Fill in the values for lonvalue_1d and latvalue_1d of rectilinear domain
457     Range of longitude value from 0 - 360
458     Range of latitude value from -90 - +90
459   */
460   void CDomain::fillInRectilinearLonLat()
461   {
462     if (!lonvalue_rectilinear_read_from_file.isEmpty() && lonvalue_2d.isEmpty() && lonvalue_1d.isEmpty())
463     {
464       lonvalue_1d.resize(ni);
465       for (int idx = 0; idx < ni; ++idx)
466         lonvalue_1d(idx) = lonvalue_rectilinear_read_from_file(idx+ibegin);
467       lon_start.setValue(lonvalue_rectilinear_read_from_file(0));
468       lon_end.setValue(lonvalue_rectilinear_read_from_file(ni_glo-1));
469     }
470     else if (!hasLonInReadFile_)
471     {
472       if (!lonvalue_2d.isEmpty()) lonvalue_2d.free();
473       lonvalue_1d.resize(ni);
474       double lonRange = lon_end - lon_start;
475       double lonStep = (1 == ni_glo.getValue()) ? lonRange : lonRange/double(ni_glo.getValue()-1);
476
477        // Assign lon value
478       for (int i = 0; i < ni; ++i)
479       {
480         if (0 == (ibegin + i))
481         {
482           lonvalue_1d(i) = lon_start;
483         }
484         else if (ni_glo == (ibegin + i + 1))
485         {
486           lonvalue_1d(i) = lon_end;
487         }
488         else
489         {
490           lonvalue_1d(i) = (ibegin + i) * lonStep  + lon_start;
491         }
492       }
493     }
494
495
496     if (!latvalue_rectilinear_read_from_file.isEmpty() && latvalue_2d.isEmpty() && latvalue_1d.isEmpty())
497     {
498       latvalue_1d.resize(nj);
499       for (int idx = 0; idx < nj; ++idx)
500         latvalue_1d(idx) = latvalue_rectilinear_read_from_file(idx+jbegin);
501       lat_start.setValue(latvalue_rectilinear_read_from_file(0));
502       lat_end.setValue(latvalue_rectilinear_read_from_file(nj_glo-1));
503     }
504     else if (!hasLatInReadFile_)
505     {
506       if (!latvalue_2d.isEmpty()) latvalue_1d.free();
507       latvalue_1d.resize(nj);
508
509       double latRange = lat_end - lat_start;
510       double latStep = (1 == nj_glo.getValue()) ? latRange : latRange/double(nj_glo.getValue()-1);
511
512       for (int j = 0; j < nj; ++j)
513       {
514         if (0 == (jbegin + j))
515         {
516            latvalue_1d(j) = lat_start;
517         }
518         else if (nj_glo == (jbegin + j + 1))
519         {
520            latvalue_1d(j) = lat_end;
521         }
522         else
523         {
524           latvalue_1d(j) =  (jbegin + j) * latStep + lat_start;
525         }
526       }
527     }
528   }
529
530    /*
531      Fill in 2D longitude and latitude of curvilinear domain read from a file.
532      If there are already longitude and latitude defined by model. We just ignore read value.
533    */
534   void CDomain::fillInCurvilinearLonLat()
535   {
536     if (!lonvalue_curvilinear_read_from_file.isEmpty() && lonvalue_2d.isEmpty() && lonvalue_1d.isEmpty())
537     {
538       lonvalue_2d.resize(ni,nj);
539       for (int jdx = 0; jdx < nj; ++jdx)
540        for (int idx = 0; idx < ni; ++idx)
541         lonvalue_2d(idx,jdx) = lonvalue_curvilinear_read_from_file(idx, jdx);
542
543       lonvalue_curvilinear_read_from_file.free();
544     }
545
546     if (!latvalue_curvilinear_read_from_file.isEmpty() && latvalue_2d.isEmpty() && latvalue_1d.isEmpty())
547     {
548       latvalue_2d.resize(ni,nj);
549       for (int jdx = 0; jdx < nj; ++jdx)
550        for (int idx = 0; idx < ni; ++idx)
551           latvalue_2d(idx,jdx) = latvalue_curvilinear_read_from_file(idx, jdx);
552
553       latvalue_curvilinear_read_from_file.free();
554     }
555
556     if (!bounds_lonvalue_curvilinear_read_from_file.isEmpty() && bounds_lon_2d.isEmpty() && bounds_lon_1d.isEmpty())
557     {
558       bounds_lon_2d.resize(nvertex,ni,nj);
559       for (int jdx = 0; jdx < nj; ++jdx)
560        for (int idx = 0; idx < ni; ++idx)
561          for (int ndx = 0; ndx < nvertex; ++ndx)
562            bounds_lon_2d(ndx,idx,jdx) = bounds_lonvalue_curvilinear_read_from_file(ndx,idx, jdx);
563
564       bounds_lonvalue_curvilinear_read_from_file.free();
565     }
566
567     if (!bounds_latvalue_curvilinear_read_from_file.isEmpty() && bounds_lat_2d.isEmpty() && bounds_lat_1d.isEmpty())
568     {
569       bounds_lat_2d.resize(nvertex,ni,nj);
570       for (int jdx = 0; jdx < nj; ++jdx)
571        for (int idx = 0; idx < ni; ++idx)
572          for (int ndx = 0; ndx < nvertex; ++ndx)
573            bounds_lat_2d(ndx,idx,jdx) = bounds_latvalue_curvilinear_read_from_file(ndx,idx, jdx);
574
575       bounds_latvalue_curvilinear_read_from_file.free();
576     }
577
578   }
579
580    /*
581      Fill in longitude and latitude of unstructured domain read from a file
582      If there are already longitude and latitude defined by model. We just igonore reading value.
583    */
584   void CDomain::fillInUnstructuredLonLat()
585   {
586     if (i_index.isEmpty())
587     {
588       i_index.resize(ni);
589       for(int idx = 0; idx < ni; ++idx) i_index(idx)=ibegin+idx;
590     }
591
592     if (!lonvalue_unstructured_read_from_file.isEmpty() && lonvalue_1d.isEmpty())
593     {
594        lonvalue_1d.resize(ni);
595        for (int idx = 0; idx < ni; ++idx)
596          lonvalue_1d(idx) = lonvalue_unstructured_read_from_file(idx);
597
598        // We dont need these values anymore, so just delete them
599        lonvalue_unstructured_read_from_file.free();
600     } 
601
602     if (!latvalue_unstructured_read_from_file.isEmpty() && latvalue_1d.isEmpty())
603     {
604        latvalue_1d.resize(ni);
605        for (int idx = 0; idx < ni; ++idx)
606          latvalue_1d(idx) =  latvalue_unstructured_read_from_file(idx);
607
608        // We dont need these values anymore, so just delete them
609        latvalue_unstructured_read_from_file.free();
610     }
611
612     if (!bounds_lonvalue_unstructured_read_from_file.isEmpty() && bounds_lon_1d.isEmpty())
613     {
614        int nbVertex = nvertex;
615        bounds_lon_1d.resize(nbVertex,ni);
616        for (int idx = 0; idx < ni; ++idx)
617          for (int jdx = 0; jdx < nbVertex; ++jdx)
618            bounds_lon_1d(jdx,idx) = bounds_lonvalue_unstructured_read_from_file(jdx, idx);
619
620        // We dont need these values anymore, so just delete them
621        bounds_lonvalue_unstructured_read_from_file.free();
622     }
623
624     if (!bounds_latvalue_unstructured_read_from_file.isEmpty() && bounds_lat_1d.isEmpty())
625     {
626        int nbVertex = nvertex;
627        bounds_lat_1d.resize(nbVertex,ni);
628        for (int idx = 0; idx < ni; ++idx)
629          for (int jdx = 0; jdx < nbVertex; ++jdx)
630            bounds_lat_1d(jdx,idx) = bounds_latvalue_unstructured_read_from_file(jdx, idx);
631
632        // We dont need these values anymore, so just delete them
633        bounds_latvalue_unstructured_read_from_file.free();
634     }
635   }
636
637  /*
638    Get global longitude and latitude of rectilinear domain.
639  */
640   void CDomain::AllgatherRectilinearLonLat(CArray<double,1>& lon, CArray<double,1>& lat, CArray<double,1>& lon_g, CArray<double,1>& lat_g)
641   {
642     CContext* context = CContext::getCurrent();
643    // For now the assumption is that secondary server pools consist of the same number of procs.
644    // CHANGE the line below if the assumption changes.
645    CContextClient* client = (0 != context->clientPrimServer.size()) ? context->clientPrimServer[0] : context->client;
646     lon_g.resize(ni_glo) ;
647     lat_g.resize(nj_glo) ;
648
649
650     int* ibegin_g = new int[client->clientSize] ;
651     int* jbegin_g = new int[client->clientSize] ;
652     int* ni_g = new int[client->clientSize] ;
653     int* nj_g = new int[client->clientSize] ;
654     int v ;
655     v=ibegin ;
656     MPI_Allgather(&v,1,MPI_INT,ibegin_g,1,MPI_INT,client->intraComm) ;
657     v=jbegin ;
658     MPI_Allgather(&v,1,MPI_INT,jbegin_g,1,MPI_INT,client->intraComm) ;
659     v=ni ;
660     MPI_Allgather(&v,1,MPI_INT,ni_g,1,MPI_INT,client->intraComm) ;
661     v=nj ;
662     MPI_Allgather(&v,1,MPI_INT,nj_g,1,MPI_INT,client->intraComm) ;
663
664     MPI_Allgatherv(lon.dataFirst(),ni,MPI_DOUBLE,lon_g.dataFirst(),ni_g, ibegin_g,MPI_DOUBLE,client->intraComm) ;
665     MPI_Allgatherv(lat.dataFirst(),nj,MPI_DOUBLE,lat_g.dataFirst(),nj_g, jbegin_g,MPI_DOUBLE,client->intraComm) ;
666
667      delete[] ibegin_g ;
668      delete[] jbegin_g ;
669      delete[] ni_g ;
670      delete[] nj_g ;
671   }
672
673   void CDomain::fillInRectilinearBoundLonLat(CArray<double,1>& lon, CArray<double,1>& lat,
674                                              CArray<double,2>& boundsLon, CArray<double,2>& boundsLat)
675   {
676     int i,j,k;
677
678     const int nvertexValue = 4;
679     boundsLon.resize(nvertexValue,ni*nj);
680
681     if (ni_glo>1)
682     {
683       double lonStepStart = lon(1)-lon(0);
684       bounds_lon_start=lon(0) - lonStepStart/2;
685       double lonStepEnd = lon(ni_glo-1)-lon(ni_glo-2);
686       bounds_lon_end=lon(ni_glo-1) + lonStepEnd/2;
687       double errorBoundsLon = std::abs(360-std::abs(bounds_lon_end-bounds_lon_start));
688
689       // if errorBoundsLon is reasonably small (0.1 x cell size) consider it as closed in longitude
690       if (errorBoundsLon < std::abs(lonStepStart)*1e-1 || errorBoundsLon < std::abs(lonStepEnd)*1e-1 )
691       {
692         bounds_lon_start= (lon(0) + lon(ni_glo-1)-360)/2 ;
693         bounds_lon_end= (lon(0) +360 + lon(ni_glo-1))/2 ;
694       }
695     }
696     else
697     {
698       if (bounds_lon_start.isEmpty()) bounds_lon_start=-180. ;
699       if (bounds_lon_end.isEmpty()) bounds_lon_end=180.-1e-8 ;
700     }
701
702     for(j=0;j<nj;++j)
703       for(i=0;i<ni;++i)
704       {
705         k=j*ni+i;
706         boundsLon(0,k) = boundsLon(1,k) = (0 == (ibegin + i)) ? bounds_lon_start
707                                                               : (lon(ibegin + i)+lon(ibegin + i-1))/2;
708         boundsLon(2,k) = boundsLon(3,k) = ((ibegin + i + 1) == ni_glo) ? bounds_lon_end
709                                                                        : (lon(ibegin + i + 1)+lon(ibegin + i))/2;
710       }
711
712
713    boundsLat.resize(nvertexValue,nj*ni);
714    bool isNorthPole=false ;
715    bool isSouthPole=false ;
716    if (std::abs(90 - std::abs(lat(0))) < NumTraits<double>::epsilon()) isNorthPole = true;
717    if (std::abs(90 - std::abs(lat(nj_glo-1))) < NumTraits<double>::epsilon()) isSouthPole = true;
718
719    // lat boundaries beyond pole the assimilate it to pole
720    // lat boundarie is relativelly close to pole (0.1 x cell size) assimilate it to pole
721    if (nj_glo>1)
722    {
723      double latStepStart = lat(1)-lat(0);
724      if (isNorthPole) bounds_lat_start=lat(0);
725      else
726      {
727        bounds_lat_start=lat(0)-latStepStart/2;
728        if (bounds_lat_start >= 90 ) bounds_lat_start=90 ;
729        else if (bounds_lat_start <= -90 ) bounds_lat_start=-90 ;
730        else if (bounds_lat_start <= 90 && bounds_lat_start >= lat(0))
731        {
732          if ( std::abs(90-bounds_lat_start) <= 0.1*std::abs(latStepStart)) bounds_lat_start=90 ;
733        }
734        else if (bounds_lat_start >= -90 && bounds_lat_start <= lat(0))
735        {
736          if ( std::abs(-90 - bounds_lat_start) <= 0.1*std::abs(latStepStart)) bounds_lat_start=-90 ;
737        }
738      }
739
740      double latStepEnd = lat(nj_glo-1)-lat(nj_glo-2);
741      if (isSouthPole) bounds_lat_end=lat(nj_glo-1);
742      else
743      {
744        bounds_lat_end=lat(nj_glo-1)+latStepEnd/2;
745
746        if (bounds_lat_end >= 90 ) bounds_lat_end=90 ;
747        else if (bounds_lat_end <= -90 ) bounds_lat_end=-90 ;
748        else if (bounds_lat_end <= 90 && bounds_lat_end >= lat(nj_glo-1))
749        {
750          if ( std::abs(90-bounds_lat_end) <= 0.1*std::abs(latStepEnd)) bounds_lat_end=90 ;
751        }
752        else if (bounds_lat_end >= -90 && bounds_lat_end <= lat(nj_glo-1))
753        {
754          if ( std::abs(-90 - bounds_lat_end) <= 0.1*std::abs(latStepEnd)) bounds_lat_end=-90 ;
755        }
756      }
757    }
758    else
759    {
760      if (bounds_lat_start.isEmpty()) bounds_lat_start=-90. ;
761      if (bounds_lat_end.isEmpty()) bounds_lat_end=90 ;
762    }
763
764    for(j=0;j<nj;++j)
765      for(i=0;i<ni;++i)
766      {
767        k=j*ni+i;
768        boundsLat(1,k) = boundsLat(2,k) = (0 == (jbegin + j)) ? bounds_lat_start
769                                                              : (lat(jbegin + j)+lat(jbegin + j-1))/2;
770        boundsLat(0,k) = boundsLat(3,k) = ((jbegin + j +1) == nj_glo) ? bounds_lat_end
771                                                                      : (lat(jbegin + j + 1)+lat(jbegin + j))/2;
772      }
773   }
774
775   /*
776     General check of the domain to verify its mandatory attributes
777   */
778   void CDomain::checkDomain(void)
779   {
780     if (type.isEmpty())
781     {
782       ERROR("CDomain::checkDomain(void)",
783             << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
784             << "The domain type is mandatory, "
785             << "please define the 'type' attribute.")
786     }
787
788     if (type == type_attr::gaussian) 
789     {
790        hasPole=true ;
791        type.setValue(type_attr::unstructured) ;
792      }
793      else if (type == type_attr::rectilinear) hasPole=true ;
794
795     if (type == type_attr::unstructured)
796     {
797        if (ni_glo.isEmpty())
798        {
799          ERROR("CDomain::checkDomain(void)",
800                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
801                << "The global domain is badly defined, "
802                << "the mandatory 'ni_glo' attribute is missing.")
803        }
804        else if (ni_glo <= 0)
805        {
806          ERROR("CDomain::checkDomain(void)",
807                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
808                << "The global domain is badly defined, "
809                << "'ni_glo' attribute should be strictly positive so 'ni_glo = " << ni_glo.getValue() << "' is invalid.")
810        }
811        isUnstructed_ = true;
812        nj_glo = 1;
813        nj = 1;
814        jbegin = 0;
815        if (!i_index.isEmpty()) ni = i_index.numElements();
816        j_index.resize(ni);
817        for(int i=0;i<ni;++i) j_index(i)=0;
818
819        if (!area.isEmpty())
820          area.transposeSelf(1, 0);
821     }
822
823     if (ni_glo.isEmpty())
824     {
825       ERROR("CDomain::checkDomain(void)",
826             << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
827             << "The global domain is badly defined, "
828             << "the mandatory 'ni_glo' attribute is missing.")
829     }
830     else if (ni_glo <= 0)
831     {
832       ERROR("CDomain::checkDomain(void)",
833             << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
834             << "The global domain is badly defined, "
835             << "'ni_glo' attribute should be strictly positive so 'ni_glo = " << ni_glo.getValue() << "' is invalid.")
836     }
837
838     if (nj_glo.isEmpty())
839     {
840       ERROR("CDomain::checkDomain(void)",
841             << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
842             << "The global domain is badly defined, "
843             << "the mandatory 'nj_glo' attribute is missing.")
844     }
845     else if (nj_glo <= 0)
846     {
847       ERROR("CDomain::checkDomain(void)",
848             << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
849             << "The global domain is badly defined, "
850             << "'nj_glo' attribute should be strictly positive so 'nj_glo = " << nj_glo.getValue() << "' is invalid.")
851     }
852
853     checkLocalIDomain();
854     checkLocalJDomain();
855
856     if (i_index.isEmpty())
857     {
858       i_index.resize(ni*nj);
859       for (int j = 0; j < nj; ++j)
860         for (int i = 0; i < ni; ++i) i_index(i+j*ni) = i+ibegin;
861     }
862
863     if (j_index.isEmpty())
864     {
865       j_index.resize(ni*nj);
866       for (int j = 0; j < nj; ++j)
867         for (int i = 0; i < ni; ++i) j_index(i+j*ni) = j+jbegin;
868     }
869   }
870
871   size_t CDomain::getGlobalWrittenSize(void)
872   {
873     return ni_glo*nj_glo ;
874   }
875   //----------------------------------------------------------------
876
877   // Check validity of local domain on using the combination of 3 parameters: ibegin, ni and i_index
878   void CDomain::checkLocalIDomain(void)
879   {
880      // If ibegin and ni are provided then we use them to check the validity of local domain
881      if (i_index.isEmpty() && !ibegin.isEmpty() && !ni.isEmpty())
882      {
883        if ((ni.getValue() < 0 || ibegin.getValue() < 0) || ((ibegin.getValue() + ni.getValue()) > ni_glo.getValue()))
884        {
885          ERROR("CDomain::checkLocalIDomain(void)",
886                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
887                << "The local domain is wrongly defined,"
888                << " check the attributes 'ni_glo' (" << ni_glo.getValue() << "), 'ni' (" << ni.getValue() << ") and 'ibegin' (" << ibegin.getValue() << ")");
889        }
890      }
891
892      // i_index has higher priority than ibegin and ni
893      if (!i_index.isEmpty())
894      {
895        int minIIndex = (0 < i_index.numElements()) ? i_index(0) : 0;
896        if (ni.isEmpty()) 
897        {         
898         // No information about ni
899          int minIndex = ni_glo - 1;
900          int maxIndex = 0;
901          for (int idx = 0; idx < i_index.numElements(); ++idx)
902          {
903            if (i_index(idx) < minIndex) minIndex = i_index(idx);
904            if (i_index(idx) > maxIndex) maxIndex = i_index(idx);
905          }
906          ni = maxIndex - minIndex + 1; 
907          minIIndex = minIIndex;         
908        }
909
910        // It's not so correct but if ibegin is not the first value of i_index
911        // then data on local domain has user-defined distribution. In this case, ibegin, ni have no meaning.
912        if (ibegin.isEmpty()) ibegin = minIIndex;
913      }
914      else if (ibegin.isEmpty() && ni.isEmpty())
915      {
916        ibegin = 0;
917        ni = ni_glo;
918      }
919      else if ((!ibegin.isEmpty() && ni.isEmpty()) || (ibegin.isEmpty() && !ni.isEmpty()))
920      {
921        ERROR("CDomain::checkLocalIDomain(void)",
922              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
923              << "The local domain is wrongly defined," << endl
924              << "i_index is empty and either 'ni' or 'ibegin' is not defined. " 
925              << "If 'ni' and 'ibegin' are used to define a domain, both of them must not be empty.");
926      }
927       
928
929      if ((ni.getValue() < 0 || ibegin.getValue() < 0))
930      {
931        ERROR("CDomain::checkLocalIDomain(void)",
932              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
933              << "The local domain is wrongly defined,"
934              << " check the attributes 'ni_glo' (" << ni_glo.getValue() << "), 'ni' (" << ni.getValue() << ") and 'ibegin' (" << ibegin.getValue() << ")");
935      }
936   }
937
938   // Check validity of local domain on using the combination of 3 parameters: jbegin, nj and j_index
939   void CDomain::checkLocalJDomain(void)
940   {
941    // If jbegin and nj are provided then we use them to check the validity of local domain
942     if (j_index.isEmpty() && !jbegin.isEmpty() && !nj.isEmpty())
943     {
944       if ((nj.getValue() < 0 || jbegin.getValue() < 0) || (jbegin.getValue() + nj.getValue()) > nj_glo.getValue())
945       {
946         ERROR("CDomain::checkLocalJDomain(void)",
947                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
948                << "The local domain is wrongly defined,"
949                << " check the attributes 'nj_glo' (" << nj_glo.getValue() << "), 'nj' (" << nj.getValue() << ") and 'jbegin' (" << jbegin.getValue() << ")");
950       }
951     }
952
953     if (!j_index.isEmpty())
954     {
955        int minJIndex = (0 < j_index.numElements()) ? j_index(0) : 0;
956        if (nj.isEmpty()) 
957        {
958          // No information about nj
959          int minIndex = nj_glo - 1;
960          int maxIndex = 0;
961          for (int idx = 0; idx < j_index.numElements(); ++idx)
962          {
963            if (j_index(idx) < minIndex) minIndex = j_index(idx);
964            if (j_index(idx) > maxIndex) maxIndex = j_index(idx);
965          }
966          nj = maxIndex - minIndex + 1;
967          minJIndex = minIndex; 
968        } 
969        // It's the same as checkLocalIDomain. It's not so correct but if jbegin is not the first value of j_index
970        // then data on local domain has user-defined distribution. In this case, jbegin has no meaning.
971       if (jbegin.isEmpty()) jbegin = minJIndex;       
972     }
973     else if (jbegin.isEmpty() && nj.isEmpty())
974     {
975       jbegin = 0;
976       nj = nj_glo;
977     }     
978
979
980     if ((nj.getValue() < 0 || jbegin.getValue() < 0))
981     {
982       ERROR("CDomain::checkLocalJDomain(void)",
983              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
984              << "The local domain is wrongly defined,"
985              << " check the attributes 'nj_glo' (" << nj_glo.getValue() << "), 'nj' (" << nj.getValue() << ") and 'jbegin' (" << jbegin.getValue() << ")");
986     }
987   }
988
989   //----------------------------------------------------------------
990
991   void CDomain::checkMask(void)
992   {
993      if (!mask_1d.isEmpty() && !mask_2d.isEmpty())
994        ERROR("CDomain::checkMask(void)",
995              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
996              << "Both mask_1d and mask_2d are defined but only one can be used at the same time." << std::endl
997              << "Please define only one mask: 'mask_1d' or 'mask_2d'.");
998
999      if (!mask_1d.isEmpty() && mask_2d.isEmpty())
1000      {
1001        if (mask_1d.numElements() != i_index.numElements())
1002          ERROR("CDomain::checkMask(void)",
1003                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1004                << "'mask_1d' does not have the same size as the local domain." << std::endl
1005                << "Local size is " << i_index.numElements() << "." << std::endl
1006                << "Mask size is " << mask_1d.numElements() << ".");
1007      }
1008
1009      if (mask_1d.isEmpty() && !mask_2d.isEmpty())
1010      {
1011        if (mask_2d.extent(0) != ni || mask_2d.extent(1) != nj)
1012          ERROR("CDomain::checkMask(void)",
1013                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1014                << "The mask does not have the same size as the local domain." << std::endl
1015                << "Local size is " << ni.getValue() << " x " << nj.getValue() << "." << std::endl
1016                << "Mask size is " << mask_2d.extent(0) << " x " << mask_2d.extent(1) << ".");
1017      }
1018
1019      if (!mask_2d.isEmpty())
1020      {
1021        domainMask.resize(mask_2d.extent(0) * mask_2d.extent(1));
1022        for (int j = 0; j < nj; ++j)
1023          for (int i = 0; i < ni; ++i) domainMask(i+j*ni) = mask_2d(i,j);
1024//        mask_2d.reset();
1025      }
1026      else if (mask_1d.isEmpty())
1027      {
1028        domainMask.resize(i_index.numElements());
1029        for (int i = 0; i < i_index.numElements(); ++i) domainMask(i) = true;
1030      }
1031      else
1032      {
1033      domainMask.resize(mask_1d.numElements());
1034      domainMask=mask_1d ;
1035     }
1036   }
1037
1038   //----------------------------------------------------------------
1039
1040   void CDomain::checkDomainData(void)
1041   {
1042      if (data_dim.isEmpty())
1043      {
1044        data_dim.setValue(1);
1045      }
1046      else if (!(data_dim.getValue() == 1 || data_dim.getValue() == 2))
1047      {
1048        ERROR("CDomain::checkDomainData(void)",
1049              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1050              << "The data dimension is invalid, 'data_dim' must be 1 or 2 not << " << data_dim.getValue() << ".");
1051      }
1052
1053      if (data_ibegin.isEmpty())
1054         data_ibegin.setValue(0);
1055      if (data_jbegin.isEmpty())
1056         data_jbegin.setValue(0);
1057
1058      if (data_ni.isEmpty())
1059      {
1060        data_ni.setValue((data_dim == 1) ? (ni.getValue() * nj.getValue()) : ni.getValue());
1061      }
1062      else if (data_ni.getValue() < 0)
1063      {
1064        ERROR("CDomain::checkDomainData(void)",
1065              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1066              << "The data size cannot be negative ('data_ni' = " << data_ni.getValue() << ").");
1067      }
1068
1069      if (data_nj.isEmpty())
1070      {
1071        data_nj.setValue((data_dim.getValue() == 1) ? (ni.getValue() * nj.getValue()) : nj.getValue());
1072      }
1073      else if (data_nj.getValue() < 0)
1074      {
1075        ERROR("CDomain::checkDomainData(void)",
1076              << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1077              << "The data size cannot be negative ('data_nj' = " << data_nj.getValue() << ").");
1078      }
1079   }
1080
1081   //----------------------------------------------------------------
1082
1083   void CDomain::checkCompression(void)
1084   {
1085      if (!data_i_index.isEmpty())
1086      {
1087        if (!data_j_index.isEmpty() &&
1088            data_j_index.numElements() != data_i_index.numElements())
1089        {
1090           ERROR("CDomain::checkCompression(void)",
1091                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1092                 << "'data_i_index' and 'data_j_index' arrays must have the same size." << std::endl
1093                 << "'data_i_index' size = " << data_i_index.numElements() << std::endl
1094                 << "'data_j_index' size = " << data_j_index.numElements());
1095        }
1096
1097        if (2 == data_dim)
1098        {
1099          if (data_j_index.isEmpty())
1100          {
1101             ERROR("CDomain::checkCompression(void)",
1102                   << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1103                   << "'data_j_index' must be defined when 'data_i_index' is set and 'data_dim' is 2.");
1104          }
1105        }
1106        else // (1 == data_dim)
1107        {
1108          if (data_j_index.isEmpty())
1109          {
1110            data_j_index.resize(data_ni);
1111            for (int j = 0; j < data_ni; ++j) data_j_index(j) = 0;
1112          }
1113        }
1114      }
1115      else
1116      {
1117        if (data_dim == 2 && !data_j_index.isEmpty())
1118          ERROR("CDomain::checkCompression(void)",
1119                << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1120                << "'data_i_index' must be defined when 'data_j_index' is set and 'data_dim' is 2.");
1121
1122        if (1 == data_dim)
1123        {
1124          data_i_index.resize(data_ni);
1125          data_j_index.resize(data_ni);
1126
1127          for (int i = 0; i < data_ni; ++i)
1128          {
1129            data_i_index(i) = i;
1130            data_j_index(i) = 0;
1131          }
1132        }
1133        else // (data_dim == 2)
1134        {
1135          const int dsize = data_ni * data_nj;
1136          data_i_index.resize(dsize);
1137          data_j_index.resize(dsize);
1138
1139          for(int count = 0, j = 0; j < data_nj; ++j)
1140          {
1141            for(int i = 0; i < data_ni; ++i, ++count)
1142            {
1143              data_i_index(count) = i;
1144              data_j_index(count) = j;
1145            }
1146          }
1147        }
1148      }
1149   }
1150
1151   //----------------------------------------------------------------
1152   void CDomain::computeLocalMask(void)
1153   {
1154     localMask.resize(i_index.numElements()) ;
1155     localMask=false ;
1156
1157     size_t dn=data_i_index.numElements() ;
1158     int i,j ;
1159     size_t k,ind ;
1160
1161     for(k=0;k<dn;k++)
1162     {
1163       if (data_dim==2)
1164       {
1165          i=data_i_index(k)+data_ibegin ;
1166          j=data_j_index(k)+data_jbegin ;
1167          if (i>=0 && i<ni && j>=0 && j<nj)
1168          {
1169            ind=j*ni+i ;
1170            localMask(ind)=domainMask(ind) ;
1171          }
1172       }
1173       else
1174       {
1175          i=data_i_index(k)+data_ibegin ;
1176          if (i>=0 && i<i_index.numElements())
1177          {
1178            ind=i ;
1179            localMask(ind)=domainMask(ind) ;
1180          }
1181       }
1182     }
1183   }
1184
1185   void CDomain::checkEligibilityForCompressedOutput(void)
1186   {
1187     // We don't check if the mask or the indexes are valid here, just if they have been defined at this point.
1188     isCompressible_ = !mask_1d.isEmpty() || !mask_2d.isEmpty() || !data_i_index.isEmpty();
1189   }
1190
1191   //----------------------------------------------------------------
1192
1193   /*
1194     Fill in longitude, latitude, bounds, and area into internal values (lonvalue, latvalue, bounds_lonvalue, bounds_latvalue, areavalue)
1195     which will be used by XIOS.
1196   */
1197   void CDomain::completeLonLatClient(void)
1198   {
1199     bool lonlatValueExisted = (0 != lonvalue.numElements()) || (0 != latvalue.numElements());
1200     checkBounds() ;
1201     checkArea() ;
1202
1203     if (!lonvalue_2d.isEmpty() && !lonlatValueExisted)
1204     {
1205       lonvalue.resize(ni * nj);
1206       latvalue.resize(ni * nj);
1207       if (hasBounds)
1208       {
1209         bounds_lonvalue.resize(nvertex, ni * nj);
1210         bounds_latvalue.resize(nvertex, ni * nj);
1211       }
1212
1213       for (int j = 0; j < nj; ++j)
1214       {
1215         for (int i = 0; i < ni; ++i)
1216         {
1217           int k = j * ni + i;
1218
1219           lonvalue(k) = lonvalue_2d(i,j);
1220           latvalue(k) = latvalue_2d(i,j);
1221
1222           if (hasBounds)
1223           {
1224             for (int n = 0; n < nvertex; ++n)
1225             {
1226               bounds_lonvalue(n,k) = bounds_lon_2d(n,i,j);
1227               bounds_latvalue(n,k) = bounds_lat_2d(n,i,j);
1228             }
1229           }
1230         }
1231       }
1232     }
1233     else if (!lonvalue_1d.isEmpty()  && !lonlatValueExisted)
1234     {
1235       if (type_attr::rectilinear == type)
1236       {
1237         if (ni == lonvalue_1d.numElements() && nj == latvalue_1d.numElements())
1238         {
1239           lonvalue.resize(ni * nj);
1240           latvalue.resize(ni * nj);
1241           if (hasBounds)
1242           {
1243             bounds_lonvalue.resize(nvertex, ni * nj);
1244             bounds_latvalue.resize(nvertex, ni * nj);
1245           }
1246
1247           for (int j = 0; j < nj; ++j)
1248           {
1249             for (int i = 0; i < ni; ++i)
1250             {
1251               int k = j * ni + i;
1252
1253               lonvalue(k) = lonvalue_1d(i);
1254               latvalue(k) = latvalue_1d(j);
1255
1256               if (hasBounds)
1257               {
1258                 for (int n = 0; n < nvertex; ++n)
1259                 {
1260                   bounds_lonvalue(n,k) = bounds_lon_1d(n,i);
1261                   bounds_latvalue(n,k) = bounds_lat_1d(n,j);
1262                 }
1263               }
1264             }
1265           }
1266         }
1267         else if (i_index.numElements() == lonvalue_1d.numElements() && j_index.numElements() == latvalue_1d.numElements()  && !lonlatValueExisted)
1268         {
1269           lonvalue.reference(lonvalue_1d);
1270           latvalue.reference(latvalue_1d);
1271            if (hasBounds)
1272           {
1273             bounds_lonvalue.reference(bounds_lon_1d);
1274             bounds_latvalue.reference(bounds_lat_1d);
1275           }
1276         }
1277         else
1278           ERROR("CDomain::completeLonClient(void)",
1279                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1280                 << "'lonvalue_1d' and 'latvalue_1d' does not have the same size as the local domain." << std::endl
1281                 << "'lonvalue_1d' size is " << lonvalue_1d.numElements() 
1282                 << " and 'latvalue_1d' size is " << latvalue_1d.numElements() << std::endl
1283                 << " They should be correspondingly " << ni.getValue() << " and "  << nj.getValue() << " or " << std::endl
1284                 << i_index.numElements() << " and "  << j_index.numElements() << ".");
1285       }
1286       else if (type == type_attr::curvilinear || type == type_attr::unstructured  && !lonlatValueExisted)
1287       {
1288         lonvalue.reference(lonvalue_1d);
1289         latvalue.reference(latvalue_1d);
1290         if (hasBounds)
1291         {
1292           bounds_lonvalue.reference(bounds_lon_1d);
1293           bounds_latvalue.reference(bounds_lat_1d);
1294         }
1295       }
1296     }
1297
1298     if (!area.isEmpty() && areavalue.isEmpty())
1299     {
1300        areavalue.resize(ni*nj);
1301       for (int j = 0; j < nj; ++j)
1302       {
1303         for (int i = 0; i < ni; ++i)
1304         {
1305           int k = j * ni + i;
1306           areavalue(k) = area(i,j);
1307         }
1308       }
1309     }
1310   }
1311
1312   /*
1313     Convert internal longitude latitude value used by XIOS to "lonvalue_*" which can be retrieved with Fortran interface
1314   */
1315   void CDomain::convertLonLatValue(void)
1316   {
1317     bool lonlatValueExisted = (0 != lonvalue.numElements()) || (0 != latvalue.numElements());
1318     if (!lonvalue_2d.isEmpty() && lonlatValueExisted)
1319     {
1320       lonvalue_2d.resize(ni,nj);
1321       latvalue_2d.resize(ni,nj);
1322       if (hasBounds)
1323       {
1324         bounds_lon_2d.resize(nvertex, ni, nj);
1325         bounds_lat_2d.resize(nvertex, ni, nj);
1326       }
1327
1328       for (int j = 0; j < nj; ++j)
1329       {
1330         for (int i = 0; i < ni; ++i)
1331         {
1332           int k = j * ni + i;
1333
1334           lonvalue_2d(i,j) = lonvalue(k);
1335           latvalue_2d(i,j) = latvalue(k);
1336
1337           if (hasBounds)
1338           {
1339             for (int n = 0; n < nvertex; ++n)
1340             {
1341               bounds_lon_2d(n,i,j) = bounds_lonvalue(n,k);
1342               bounds_lat_2d(n,i,j) = bounds_latvalue(n,k);
1343             }
1344           }
1345         }
1346       }
1347     }
1348     else if (!lonvalue_1d.isEmpty()  && lonlatValueExisted)
1349     {
1350       if (type_attr::rectilinear == type)
1351       {
1352         if (ni == lonvalue_1d.numElements() && nj == latvalue_1d.numElements())
1353         {
1354           lonvalue.resize(ni * nj);
1355           latvalue.resize(ni * nj);
1356           if (hasBounds)
1357           {
1358             bounds_lonvalue.resize(nvertex, ni * nj);
1359             bounds_latvalue.resize(nvertex, ni * nj);
1360           }
1361
1362           for (int j = 0; j < nj; ++j)
1363           {
1364             for (int i = 0; i < ni; ++i)
1365             {
1366               int k = j * ni + i;
1367
1368               lonvalue(k) = lonvalue_1d(i);
1369               latvalue(k) = latvalue_1d(j);
1370
1371               if (hasBounds)
1372               {
1373                 for (int n = 0; n < nvertex; ++n)
1374                 {
1375                   bounds_lonvalue(n,k) = bounds_lon_1d(n,i);
1376                   bounds_latvalue(n,k) = bounds_lat_1d(n,j);
1377                 }
1378               }
1379             }
1380           }
1381         }
1382         else if (i_index.numElements() == lonvalue_1d.numElements() && j_index.numElements() == latvalue_1d.numElements()  && !lonlatValueExisted)
1383         {
1384           lonvalue.reference(lonvalue_1d);
1385           latvalue.reference(latvalue_1d);
1386            if (hasBounds)
1387           {
1388             bounds_lonvalue.reference(bounds_lon_1d);
1389             bounds_latvalue.reference(bounds_lat_1d);
1390           }
1391         }
1392         else
1393           ERROR("CDomain::completeLonClient(void)",
1394                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1395                 << "'lonvalue_1d' and 'latvalue_1d' does not have the same size as the local domain." << std::endl
1396                 << "'lonvalue_1d' size is " << lonvalue_1d.numElements() 
1397                 << " and 'latvalue_1d' size is " << latvalue_1d.numElements() << std::endl
1398                 << " They should be correspondingly " << ni.getValue() << " and "  << nj.getValue() << " or " << std::endl
1399                 << i_index.numElements() << " and "  << j_index.numElements() << ".");
1400       }
1401       else if (type == type_attr::curvilinear || type == type_attr::unstructured  && !lonlatValueExisted)
1402       {
1403         lonvalue.reference(lonvalue_1d);
1404         latvalue.reference(latvalue_1d);
1405         if (hasBounds)
1406         {
1407           bounds_lonvalue.reference(bounds_lon_1d);
1408           bounds_latvalue.reference(bounds_lat_1d);
1409         }
1410       }
1411     }
1412   }
1413
1414
1415   void CDomain::checkBounds(void)
1416   {
1417     bool hasBoundValues = (0 != bounds_lonvalue.numElements()) || (0 != bounds_latvalue.numElements());
1418     if (!nvertex.isEmpty() && nvertex > 0 && !hasBoundValues)
1419     {
1420       if (!bounds_lon_1d.isEmpty() && !bounds_lon_2d.isEmpty())
1421         ERROR("CDomain::checkBounds(void)",
1422               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1423               << "Only one longitude boundary attribute can be used but both 'bounds_lon_1d' and 'bounds_lon_2d' are defined." << std::endl
1424               << "Define only one longitude boundary attribute: 'bounds_lon_1d' or 'bounds_lon_2d'.");
1425
1426       if (!bounds_lat_1d.isEmpty() && !bounds_lat_2d.isEmpty())
1427         ERROR("CDomain::checkBounds(void)",
1428               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1429               << "Only one latitude boundary attribute can be used but both 'bounds_lat_1d' and 'bounds_lat_2d' are defined." << std::endl
1430               << "Define only one latitude boundary attribute: 'bounds_lat_1d' or 'bounds_lat_2d'.");
1431
1432       if ((!bounds_lon_1d.isEmpty() && bounds_lat_1d.isEmpty()) || (bounds_lon_1d.isEmpty() && !bounds_lat_1d.isEmpty()))
1433       {
1434         ERROR("CDomain::checkBounds(void)",
1435               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1436               << "Only 'bounds_lon_1d' or 'bounds_lat_1d' is defined." << std::endl
1437               << "Please define either both attributes or none.");
1438       }
1439
1440       if ((!bounds_lon_2d.isEmpty() && bounds_lat_2d.isEmpty()) || (bounds_lon_2d.isEmpty() && !bounds_lat_2d.isEmpty()))
1441       {
1442         ERROR("CDomain::checkBounds(void)",
1443               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1444               << "Only 'bounds_lon_2d' or 'bounds_lat_2d' is defined." << std::endl
1445               << "Please define either both attributes or none.");
1446       }
1447
1448       if (!bounds_lon_1d.isEmpty() && nvertex.getValue() != bounds_lon_1d.extent(0))
1449         ERROR("CDomain::checkBounds(void)",
1450               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1451               << "'bounds_lon_1d' dimension is not compatible with 'nvertex'." << std::endl
1452               << "'bounds_lon_1d' dimension is " << bounds_lon_1d.extent(0)
1453               << " but nvertex is " << nvertex.getValue() << ".");
1454
1455       if (!bounds_lon_2d.isEmpty() && nvertex.getValue() != bounds_lon_2d.extent(0))
1456         ERROR("CDomain::checkBounds(void)",
1457               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1458               << "'bounds_lon_2d' dimension is not compatible with 'nvertex'." << std::endl
1459               << "'bounds_lon_2d' dimension is " << bounds_lon_2d.extent(0)
1460               << " but nvertex is " << nvertex.getValue() << ".");
1461
1462       if (!bounds_lon_1d.isEmpty() && lonvalue_1d.isEmpty())
1463         ERROR("CDomain::checkBounds(void)",
1464               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1465               << "Since 'bounds_lon_1d' is defined, 'lonvalue_1d' must be defined too." << std::endl);
1466
1467       if (!bounds_lon_2d.isEmpty() && lonvalue_2d.isEmpty())
1468         ERROR("CDomain::checkBounds(void)",
1469               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1470               << "Since 'bounds_lon_2d' is defined, 'lonvalue_2d' must be defined too." << std::endl);
1471
1472       if (!bounds_lat_1d.isEmpty() && nvertex.getValue() != bounds_lat_1d.extent(0))
1473         ERROR("CDomain::checkBounds(void)",
1474               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1475               << "'bounds_lat_1d' dimension is not compatible with 'nvertex'." << std::endl
1476               << "'bounds_lat_1d' dimension is " << bounds_lat_1d.extent(0)
1477               << " but nvertex is " << nvertex.getValue() << ".");
1478
1479       if (!bounds_lat_2d.isEmpty() && nvertex.getValue() != bounds_lat_2d.extent(0))
1480         ERROR("CDomain::checkBounds(void)",
1481               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1482               << "'bounds_lat_2d' dimension is not compatible with 'nvertex'." << std::endl
1483               << "'bounds_lat_2d' dimension is " << bounds_lat_2d.extent(0)
1484               << " but nvertex is " << nvertex.getValue() << ".");
1485
1486       if (!bounds_lat_1d.isEmpty() && latvalue_1d.isEmpty())
1487         ERROR("CDomain::checkBounds(void)",
1488               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1489               << "Since 'bounds_lat_1d' is defined, 'latvalue_1d' must be defined too." << std::endl);
1490
1491       if (!bounds_lat_2d.isEmpty() && latvalue_2d.isEmpty())
1492         ERROR("CDomain::checkBounds(void)",
1493               << "Since 'bounds_lat_2d' is defined, 'latvalue_2d' must be defined too." << std::endl);
1494
1495       // In case of reading UGRID bounds values are not required
1496       hasBounds = (!bounds_lat_1d.isEmpty() || !bounds_lat_2d.isEmpty() );
1497     }
1498     else if (hasBoundValues)
1499     {
1500       hasBounds = true;       
1501     }
1502     else
1503     {
1504       hasBounds = false;
1505     }
1506   }
1507
1508   void CDomain::checkArea(void)
1509   {
1510     bool hasAreaValue = (!areavalue.isEmpty() && 0 != areavalue.numElements());
1511     hasArea = !area.isEmpty();
1512     if (hasArea && !hasAreaValue)
1513     {
1514       if (area.extent(0) != ni || area.extent(1) != nj)
1515       {
1516         ERROR("CDomain::checkArea(void)",
1517               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1518               << "The area does not have the same size as the local domain." << std::endl
1519               << "Local size is " << ni.getValue() << " x " << nj.getValue() << "." << std::endl
1520               << "Area size is " << area.extent(0) << " x " << area.extent(1) << ".");
1521       }
1522//       if (areavalue.isEmpty())
1523//       {
1524//          areavalue.resize(ni*nj);
1525//         for (int j = 0; j < nj; ++j)
1526//         {
1527//           for (int i = 0; i < ni; ++i)
1528//           {
1529//             int k = j * ni + i;
1530//             areavalue(k) = area(i,j);
1531//           }
1532//         }
1533//       }
1534     }
1535   }
1536
1537   void CDomain::checkLonLat()
1538   {
1539     if (!hasLonLat) hasLonLat = (!latvalue_1d.isEmpty() && !lonvalue_1d.isEmpty()) ||
1540                                 (!latvalue_2d.isEmpty() && !lonvalue_2d.isEmpty());
1541     bool hasLonLatValue = (0 != lonvalue.numElements()) || (0 != latvalue.numElements());
1542     if (hasLonLat && !hasLonLatValue)
1543     {
1544       if (!lonvalue_1d.isEmpty() && !lonvalue_2d.isEmpty())
1545         ERROR("CDomain::checkLonLat()",
1546               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1547               << "Only one longitude attribute can be used but both 'lonvalue_1d' and 'lonvalue_2d' are defined." << std::endl
1548               << "Define only one longitude attribute: 'lonvalue_1d' or 'lonvalue_2d'.");
1549
1550       if (!lonvalue_1d.isEmpty() && lonvalue_2d.isEmpty())
1551       {
1552         if ((type_attr::rectilinear != type) && (lonvalue_1d.numElements() != i_index.numElements()))
1553           ERROR("CDomain::checkLonLat()",
1554                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1555                 << "'lonvalue_1d' does not have the same size as the local domain." << std::endl
1556                 << "Local size is " << i_index.numElements() << "." << std::endl
1557                 << "'lonvalue_1d' size is " << lonvalue_1d.numElements() << ".");
1558       }
1559
1560       if (lonvalue_1d.isEmpty() && !lonvalue_2d.isEmpty())
1561       {
1562         if (lonvalue_2d.extent(0) != ni || lonvalue_2d.extent(1) != nj)
1563           ERROR("CDomain::checkLonLat()",
1564                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1565                 << "'lonvalue_2d' does not have the same size as the local domain." << std::endl
1566                 << "Local size is " << ni.getValue() << " x " << nj.getValue() << "." << std::endl
1567                 << "'lonvalue_2d' size is " << lonvalue_2d.extent(0) << " x " << lonvalue_2d.extent(1) << ".");
1568       }
1569
1570       if (!latvalue_1d.isEmpty() && !latvalue_2d.isEmpty())
1571         ERROR("CDomain::checkLonLat()",
1572               << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1573               << "Only one latitude attribute can be used but both 'latvalue_1d' and 'latvalue_2d' are defined." << std::endl
1574               << "Define only one latitude attribute: 'latvalue_1d' or 'latvalue_2d'.");
1575
1576       if (!latvalue_1d.isEmpty() && latvalue_2d.isEmpty())
1577       {
1578         if ((type_attr::rectilinear != type) && (latvalue_1d.numElements() != i_index.numElements()))
1579           ERROR("CDomain::checkLonLat()",
1580                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1581                 << "'latvalue_1d' does not have the same size as the local domain." << std::endl
1582                 << "Local size is " << i_index.numElements() << "." << std::endl
1583                 << "'latvalue_1d' size is " << latvalue_1d.numElements() << ".");
1584       }
1585
1586       if (latvalue_1d.isEmpty() && !latvalue_2d.isEmpty())
1587       {
1588         if (latvalue_2d.extent(0) != ni || latvalue_2d.extent(1) != nj)
1589           ERROR("CDomain::checkLonLat()",
1590                 << "[ id = " << this->getId() << " , context = '" << CObjectFactory::GetCurrentContextId() << " ] "
1591                 << "'latvalue_2d' does not have the same size as the local domain." << std::endl
1592                 << "Local size is " << ni.getValue() << " x " << nj.getValue() << "." << std::endl
1593                 << "'latvalue_2d' size is " << latvalue_2d.extent(0) << " x " << latvalue_2d.extent(1) << ".");
1594       }
1595     }
1596   }
1597
1598   void CDomain::checkAttributesOnClientAfterTransformation()
1599   {
1600     CContext* context=CContext::getCurrent() ;
1601
1602     if (this->isClientAfterTransformationChecked) return;
1603     if (context->hasClient)
1604     {
1605      this->computeConnectedClients();
1606       if (hasLonLat)
1607         if (!context->hasServer)
1608           this->completeLonLatClient();
1609     }
1610
1611     this->isClientAfterTransformationChecked = true;
1612   }
1613
1614   //----------------------------------------------------------------
1615   // Divide function checkAttributes into 2 seperate ones
1616   // This function only checks all attributes of current domain
1617   void CDomain::checkAttributesOnClient()
1618   {
1619     if (this->isClientChecked) return;
1620     CContext* context=CContext::getCurrent();
1621
1622      if (context->hasClient && !context->hasServer)
1623      {
1624        this->checkDomain();
1625        this->checkBounds();
1626        this->checkArea();
1627        this->checkLonLat();
1628      }
1629
1630      if (context->hasClient && !context->hasServer)
1631      { // Ct client uniquement
1632         this->checkMask();
1633         this->checkDomainData();
1634         this->checkCompression();
1635         this->computeLocalMask() ;
1636      }
1637      else
1638      { // Ct serveur uniquement
1639      }
1640
1641      this->isClientChecked = true;
1642   }
1643
1644   // Send all checked attributes to server
1645   void CDomain::sendCheckedAttributes()
1646   {
1647     if (!this->isClientChecked) checkAttributesOnClient();
1648     if (!this->isClientAfterTransformationChecked) checkAttributesOnClientAfterTransformation();
1649     CContext* context=CContext::getCurrent() ;
1650
1651     if (this->isChecked) return;
1652     if (context->hasClient)
1653     {
1654       sendAttributes();
1655     }
1656     this->isChecked = true;
1657   }
1658
1659   void CDomain::checkAttributes(void)
1660   {
1661      if (this->isChecked) return;
1662      CContext* context=CContext::getCurrent() ;
1663
1664      this->checkDomain();
1665      this->checkLonLat();
1666      this->checkBounds();
1667      this->checkArea();
1668
1669      if (context->hasClient)
1670      { // Ct client uniquement
1671         this->checkMask();
1672         this->checkDomainData();
1673         this->checkCompression();
1674         this->computeLocalMask() ;
1675
1676      }
1677      else
1678      { // Ct serveur uniquement
1679      }
1680
1681      if (context->hasClient)
1682      {
1683        this->computeConnectedClients();
1684        this->completeLonLatClient();
1685      }
1686
1687      this->isChecked = true;
1688   }
1689
1690  /*!
1691     Compute the connection of a client to other clients to determine which clients to send attributes to.
1692     The sending clients are supposed to already know the distribution of receiving clients (In simple cases, it's band)
1693     The connection among clients is calculated by using global index.
1694     A client connects to other clients which holds the same global index as it.     
1695  */
1696  void CDomain::computeConnectedClients()
1697  {
1698    CContext* context=CContext::getCurrent() ;
1699   
1700    // This line should be changed soon.
1701    int nbSrvPools = (context->hasServer) ? (context->hasClient ? context->clientPrimServer.size() : 0) : 1;
1702
1703    nbSenders.clear();
1704    connectedServerRank_.clear();
1705
1706    for (int p = 0; p < nbSrvPools; ++p)
1707    {
1708      CContextClient* client = (0 != context->clientPrimServer.size()) ? context->clientPrimServer[p] : context->client;
1709      int nbServer = client->serverSize;
1710      int nbClient = client->clientSize;
1711      int rank     = client->clientRank;
1712      bool doComputeGlobalIndexServer = true;
1713
1714      if (connectedServerRank_.find(nbServer) == connectedServerRank_.end())
1715      {
1716
1717        if (indSrv_.find(nbServer) == indSrv_.end())
1718        {
1719          int i,j,i_ind,j_ind, nbIndex=i_index.numElements();
1720          int globalIndexCount = i_index.numElements();
1721          // Fill in index
1722          CArray<size_t,1> globalIndexDomain(nbIndex);
1723          size_t globalIndex;
1724
1725          for (i = 0; i < nbIndex; ++i)
1726          {
1727            i_ind=i_index(i);
1728            j_ind=j_index(i);
1729            globalIndex = i_ind + j_ind * ni_glo;
1730            globalIndexDomain(i) = globalIndex;
1731          }
1732
1733          if (globalLocalIndexMap_.empty())
1734          {
1735            for (i = 0; i < nbIndex; ++i)
1736              globalLocalIndexMap_[globalIndexDomain(i)] = i;
1737          }
1738
1739          size_t globalSizeIndex = 1, indexBegin, indexEnd;
1740          int range, clientSize = client->clientSize;
1741          std::vector<int> nGlobDomain(2);
1742          nGlobDomain[0] = this->ni_glo;
1743          nGlobDomain[1] = this->nj_glo;
1744          for (int i = 0; i < nGlobDomain.size(); ++i) globalSizeIndex *= nGlobDomain[i];
1745          indexBegin = 0;
1746          if (globalSizeIndex <= clientSize)
1747          {
1748            indexBegin = rank%globalSizeIndex;
1749            indexEnd = indexBegin;
1750          }
1751          else
1752          {
1753            for (int i = 0; i < clientSize; ++i)
1754            {
1755              range = globalSizeIndex / clientSize;
1756              if (i < (globalSizeIndex%clientSize)) ++range;
1757              if (i == client->clientRank) break;
1758              indexBegin += range;
1759            }
1760            indexEnd = indexBegin + range - 1;
1761          }
1762
1763          // Even if servers have no index, they must received something from client
1764          // We only use several client to send "empty" message to these servers
1765          CServerDistributionDescription serverDescription(nGlobDomain, nbServer);
1766          std::vector<int> serverZeroIndex;
1767          if (isUnstructed_) serverZeroIndex = serverDescription.computeServerGlobalIndexInRange(std::make_pair<size_t&,size_t&>(indexBegin, indexEnd), 0);
1768          else serverZeroIndex = serverDescription.computeServerGlobalIndexInRange(std::make_pair<size_t&,size_t&>(indexBegin, indexEnd), 1);
1769
1770          std::list<int> serverZeroIndexLeader;
1771          std::list<int> serverZeroIndexNotLeader;
1772          CContextClient::computeLeader(client->clientRank, client->clientSize, serverZeroIndex.size(), serverZeroIndexLeader, serverZeroIndexNotLeader);
1773          for (std::list<int>::iterator it = serverZeroIndexLeader.begin(); it != serverZeroIndexLeader.end(); ++it)
1774            *it = serverZeroIndex[*it];
1775
1776          CClientServerMapping* clientServerMap = new CClientServerMappingDistributed(serverDescription.getGlobalIndexRange(), client->intraComm);
1777          clientServerMap->computeServerIndexMapping(globalIndexDomain, nbServer);
1778          CClientServerMapping::GlobalIndexMap& globalIndexDomainOnServer = clientServerMap->getGlobalIndexOnServer();
1779
1780          CClientServerMapping::GlobalIndexMap::const_iterator it  = globalIndexDomainOnServer.begin(),
1781                 ite = globalIndexDomainOnServer.end();
1782          indSrv_[nbServer].swap(globalIndexDomainOnServer);
1783          connectedServerRank_[nbServer].clear();
1784          for (it = indSrv_[nbServer].begin(); it != ite; ++it)
1785            connectedServerRank_[nbServer].push_back(it->first);
1786
1787          for (std::list<int>::const_iterator it = serverZeroIndexLeader.begin(); it != serverZeroIndexLeader.end(); ++it)
1788            connectedServerRank_[nbServer].push_back(*it);
1789
1790          // Even if a client has no index, it must connect to at least one server and
1791          // send an "empty" data to this server
1792          if (connectedServerRank_[nbServer].empty())
1793            connectedServerRank_[nbServer].push_back(client->clientRank % client->serverSize);
1794
1795          nbSenders[nbServer] = clientServerMap->computeConnectedClients(client->serverSize, client->clientSize, client->intraComm, connectedServerRank_[nbServer]);
1796          delete clientServerMap;
1797        }
1798      }
1799    }
1800  }
1801
1802   /*!
1803     Compute index to write data. We only write data on the zoomed region, therefore, there should
1804     be a map between the complete grid and the reduced grid where we write data.
1805     By using global index we can easily create this kind of mapping.
1806   */
1807   void CDomain::computeWrittenIndex()
1808   { 
1809      if (computedWrittenIndex_) return;
1810      computedWrittenIndex_ = true;
1811
1812      CContext* context=CContext::getCurrent();     
1813      CContextServer* server = context->server; 
1814
1815      std::vector<int> nBegin(2), nSize(2), nBeginGlobal(2), nGlob(2);
1816      nBegin[0]       = ibegin;  nBegin[1] = jbegin;
1817      nSize[0]        = ni;      nSize[1]  = nj;
1818      nBeginGlobal[0] = 0; nBeginGlobal[1] = 0;
1819      nGlob[0]        = ni_glo;   nGlob[1] = nj_glo;
1820      CDistributionServer srvDist(server->intraCommSize, nBegin, nSize, nBeginGlobal, nGlob); 
1821      const CArray<size_t,1>& writtenGlobalIndex  = srvDist.getGlobalIndex();
1822
1823      size_t nbWritten = 0, indGlo;     
1824      std::unordered_map<size_t,size_t>::const_iterator itb = globalLocalIndexMap_.begin(),
1825                                                          ite = globalLocalIndexMap_.end(), it;         
1826      CArray<size_t,1>::const_iterator itSrvb = writtenGlobalIndex.begin(),
1827                                       itSrve = writtenGlobalIndex.end(), itSrv;
1828
1829//      for (itSrv = itSrvb; itSrv != itSrve; ++itSrv)
1830//      {
1831//        indGlo = *itSrv;
1832//        if (ite != globalLocalIndexMap_.find(indGlo))
1833//        {
1834//          ++nbWritten;
1835//        }
1836//      }
1837
1838//      localIndexToWriteOnServer.resize(nbWritten);
1839      localIndexToWriteOnServer.resize(writtenGlobalIndex.numElements());
1840
1841      nbWritten = 0;
1842      for (itSrv = itSrvb; itSrv != itSrve; ++itSrv)
1843      {
1844        indGlo = *itSrv;
1845        if (ite != globalLocalIndexMap_.find(indGlo))
1846        {
1847          localIndexToWriteOnServer(nbWritten) = globalLocalIndexMap_[indGlo];
1848          ++nbWritten;
1849        }
1850        else
1851        {
1852          localIndexToWriteOnServer(nbWritten) = 0;
1853          ++nbWritten;
1854        }
1855      }
1856     
1857      // if (isCompressible())
1858      // {
1859      //   nbWritten = 0;
1860      //   std::unordered_map<size_t,size_t> localGlobalIndexMap;
1861      //   for (itSrv = itSrvb; itSrv != itSrve; ++itSrv)
1862      //   {
1863      //     indGlo = *itSrv;
1864      //     if (ite != globalLocalIndexMap_.find(indGlo))
1865      //     {
1866      //       localGlobalIndexMap[localIndexToWriteOnServer(nbWritten)] = indGlo;
1867      //       ++nbWritten;
1868      //     }                 
1869      //   }
1870
1871      //   nbWritten = 0;
1872      //   for (int idx = 0; idx < data_i_index.numElements(); ++idx)
1873      //   {
1874      //     if (localGlobalIndexMap.end() != localGlobalIndexMap.find(data_i_index(idx)))
1875      //     {
1876      //       ++nbWritten;
1877      //     }
1878      //   }
1879
1880      //   compressedIndexToWriteOnServer.resize(nbWritten);
1881      //   nbWritten = 0;
1882      //   for (int idx = 0; idx < data_i_index.numElements(); ++idx)
1883      //   {
1884      //     if (localGlobalIndexMap.end() != localGlobalIndexMap.find(data_i_index(idx)))
1885      //     {
1886      //       compressedIndexToWriteOnServer(nbWritten) = localGlobalIndexMap[data_i_index(idx)];
1887      //       ++nbWritten;
1888      //     }
1889      //   }
1890
1891      //   numberWrittenIndexes_ = nbWritten;
1892      //   if (isDistributed())
1893      //   {           
1894      //     MPI_Allreduce(&numberWrittenIndexes_, &totalNumberWrittenIndexes_, 1, MPI_INT, MPI_SUM, server->intraComm);
1895      //     MPI_Scan(&numberWrittenIndexes_, &offsetWrittenIndexes_, 1, MPI_INT, MPI_SUM, server->intraComm);
1896      //     offsetWrittenIndexes_ -= numberWrittenIndexes_;
1897      //   }
1898      //   else
1899      //     totalNumberWrittenIndexes_ = numberWrittenIndexes_;
1900      // }     
1901   }
1902
1903  void CDomain::computeWrittenCompressedIndex(MPI_Comm writtenComm)
1904  {
1905    int writtenCommSize;
1906    MPI_Comm_size(writtenComm, &writtenCommSize);
1907    if (compressedIndexToWriteOnServer.find(writtenCommSize) != compressedIndexToWriteOnServer.end())
1908      return;
1909
1910    if (isCompressible())
1911    {
1912      size_t nbWritten = 0, indGlo;
1913      CContext* context=CContext::getCurrent();     
1914      CContextServer* server = context->server; 
1915
1916      std::vector<int> nBegin(2), nSize(2), nBeginGlobal(2), nGlob(2);
1917      nBegin[0]       = ibegin;  nBegin[1] = jbegin;
1918      nSize[0]        = ni;      nSize[1]  = nj;
1919      nBeginGlobal[0] = 0; nBeginGlobal[1] = 0;
1920      nGlob[0]        = ni_glo;   nGlob[1] = nj_glo;
1921      CDistributionServer srvDist(server->intraCommSize, nBegin, nSize, nBeginGlobal, nGlob); 
1922      const CArray<size_t,1>& writtenGlobalIndex  = srvDist.getGlobalIndex();
1923
1924      std::unordered_map<size_t,size_t>::const_iterator itb = globalLocalIndexMap_.begin(),
1925                                                          ite = globalLocalIndexMap_.end(), it;   
1926      CArray<size_t,1>::const_iterator itSrvb = writtenGlobalIndex.begin(),
1927                                       itSrve = writtenGlobalIndex.end(), itSrv;
1928      std::unordered_map<size_t,size_t> localGlobalIndexMap;
1929      for (itSrv = itSrvb; itSrv != itSrve; ++itSrv)
1930      {
1931        indGlo = *itSrv;
1932        if (ite != globalLocalIndexMap_.find(indGlo))
1933        {
1934          localGlobalIndexMap[localIndexToWriteOnServer(nbWritten)] = indGlo;
1935          ++nbWritten;
1936        }                 
1937      }
1938
1939      nbWritten = 0;
1940      for (int idx = 0; idx < data_i_index.numElements(); ++idx)
1941      {
1942        if (localGlobalIndexMap.end() != localGlobalIndexMap.find(data_i_index(idx)))
1943        {
1944          ++nbWritten;
1945        }
1946      }
1947
1948      compressedIndexToWriteOnServer[writtenCommSize].resize(nbWritten);
1949      nbWritten = 0;
1950      for (int idx = 0; idx < data_i_index.numElements(); ++idx)
1951      {
1952        if (localGlobalIndexMap.end() != localGlobalIndexMap.find(data_i_index(idx)))
1953        {
1954          compressedIndexToWriteOnServer[writtenCommSize](nbWritten) = localGlobalIndexMap[data_i_index(idx)];
1955          ++nbWritten;
1956        }
1957      }
1958
1959      numberWrittenIndexes_[writtenCommSize] = nbWritten;
1960      if (isDistributed())
1961      {
1962             
1963        MPI_Allreduce(&numberWrittenIndexes_[writtenCommSize], &totalNumberWrittenIndexes_[writtenCommSize], 1, MPI_INT, MPI_SUM, writtenComm);
1964        MPI_Scan(&numberWrittenIndexes_[writtenCommSize], &offsetWrittenIndexes_[writtenCommSize], 1, MPI_INT, MPI_SUM, writtenComm);
1965        offsetWrittenIndexes_[writtenCommSize] -= numberWrittenIndexes_[writtenCommSize];
1966      }
1967      else
1968        totalNumberWrittenIndexes_[writtenCommSize] = numberWrittenIndexes_[writtenCommSize];
1969      }
1970  }
1971
1972  /*!
1973    Send all attributes from client to connected clients
1974    The attributes will be rebuilt on receiving side
1975  */
1976  void CDomain::sendAttributes()
1977  {
1978    sendDistributionAttributes();
1979    sendIndex();       
1980    sendMask();
1981    sendLonLat();
1982    sendArea();   
1983    sendDataIndex();
1984  }
1985
1986  /*!
1987    Send global index from client to connected client(s)
1988  */
1989  void CDomain::sendIndex()
1990  {
1991    int ns, n, i, j, ind, nv, idx;
1992    std::list<CContextClient*>::iterator it;
1993    for (it=clients.begin(); it!=clients.end(); ++it)
1994    {
1995      CContextClient* client = *it;
1996
1997      int serverSize = client->serverSize;
1998      CEventClient eventIndex(getType(), EVENT_ID_INDEX);
1999
2000      list<CMessage> list_msgsIndex;
2001      list<CArray<int,1> > list_indGlob;
2002
2003      std::unordered_map<int, vector<size_t> >::const_iterator itIndex, iteIndex;
2004      iteIndex = indSrv_[serverSize].end();
2005      for (int k = 0; k < connectedServerRank_[serverSize].size(); ++k)
2006      {
2007        int nbIndGlob = 0;
2008        int rank = connectedServerRank_[serverSize][k];
2009        itIndex = indSrv_[serverSize].find(rank);
2010        if (iteIndex != itIndex)
2011          nbIndGlob = itIndex->second.size();
2012
2013        list_indGlob.push_back(CArray<int,1>(nbIndGlob));       
2014
2015        CArray<int,1>& indGlob = list_indGlob.back();
2016        for (n = 0; n < nbIndGlob; ++n)
2017        {
2018          indGlob(n) = static_cast<int>(itIndex->second[n]);
2019        }
2020
2021        list_msgsIndex.push_back(CMessage());
2022        list_msgsIndex.back() << this->getId() << (int)type; // enum ne fonctionne pour les message => ToFix
2023        list_msgsIndex.back() << isCurvilinear;
2024        list_msgsIndex.back() << list_indGlob.back(); //list_indi.back() << list_indj.back();
2025       
2026        eventIndex.push(rank, nbSenders[serverSize][rank], list_msgsIndex.back());
2027      }
2028
2029      client->sendEvent(eventIndex);
2030    }
2031  }
2032
2033  /*!
2034    Send distribution from client to other clients
2035    Because a client in a level knows correctly the grid distribution of client on the next level
2036    it calculates this distribution then sends it to the corresponding clients on the next level
2037  */
2038  void CDomain::sendDistributionAttributes(void)
2039  {
2040    std::list<CContextClient*>::iterator it;
2041    for (it=clients.begin(); it!=clients.end(); ++it)
2042    {
2043      CContextClient* client = *it;
2044      int nbServer = client->serverSize;
2045      std::vector<int> nGlobDomain(2);
2046      nGlobDomain[0] = this->ni_glo;
2047      nGlobDomain[1] = this->nj_glo;
2048
2049      CServerDistributionDescription serverDescription(nGlobDomain, nbServer);
2050      if (isUnstructed_) serverDescription.computeServerDistribution(false, 0);
2051      else serverDescription.computeServerDistribution(false, 1);
2052
2053      std::vector<std::vector<int> > serverIndexBegin = serverDescription.getServerIndexBegin();
2054      std::vector<std::vector<int> > serverDimensionSizes = serverDescription.getServerDimensionSizes();
2055
2056      CEventClient event(getType(),EVENT_ID_SERVER_ATTRIBUT);
2057      if (client->isServerLeader())
2058      {
2059        std::list<CMessage> msgs;
2060
2061        const std::list<int>& ranks = client->getRanksServerLeader();
2062        for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2063        {
2064          // Use const int to ensure CMessage holds a copy of the value instead of just a reference
2065          const int ibegin_srv = serverIndexBegin[*itRank][0];
2066          const int jbegin_srv = serverIndexBegin[*itRank][1];
2067          const int ni_srv = serverDimensionSizes[*itRank][0];
2068          const int nj_srv = serverDimensionSizes[*itRank][1];
2069
2070          msgs.push_back(CMessage());
2071          CMessage& msg = msgs.back();
2072          msg << this->getId() ;
2073          msg << isUnstructed_;
2074          msg << ni_srv << ibegin_srv << nj_srv << jbegin_srv;
2075          msg << ni_glo.getValue() << nj_glo.getValue();
2076          msg << isCompressible_;
2077
2078          event.push(*itRank,1,msg);
2079        }
2080        client->sendEvent(event);
2081      }
2082      else client->sendEvent(event);
2083    }
2084  }
2085
2086  /*!
2087    Send mask index from client to connected(s) clients   
2088  */
2089  void CDomain::sendMask()
2090  {
2091    int ns, n, i, j, ind, nv, idx;
2092    std::list<CContextClient*>::iterator it;
2093    for (it=clients.begin(); it!=clients.end(); ++it)
2094    {
2095      CContextClient* client = *it;
2096      int serverSize = client->serverSize;
2097
2098      // send area for each connected server
2099      CEventClient eventMask(getType(), EVENT_ID_MASK);
2100
2101      list<CMessage> list_msgsMask;
2102      list<CArray<bool,1> > list_mask;
2103
2104      std::unordered_map<int, vector<size_t> >::const_iterator it, iteMap;
2105      iteMap = indSrv_[serverSize].end();
2106      for (int k = 0; k < connectedServerRank_[serverSize].size(); ++k)
2107      {
2108        int nbData = 0;
2109        int rank = connectedServerRank_[serverSize][k];
2110        it = indSrv_[serverSize].find(rank);
2111        if (iteMap != it)
2112          nbData = it->second.size();
2113        list_mask.push_back(CArray<bool,1>(nbData));
2114
2115        const std::vector<size_t>& temp = it->second;
2116        for (n = 0; n < nbData; ++n)
2117        {
2118          idx = static_cast<int>(it->second[n]);
2119          list_mask.back()(n) = domainMask(globalLocalIndexMap_[idx]);
2120        }
2121
2122        list_msgsMask.push_back(CMessage());
2123        list_msgsMask.back() << this->getId() << list_mask.back();
2124        eventMask.push(rank, nbSenders[serverSize][rank], list_msgsMask.back());
2125      }
2126      client->sendEvent(eventMask);
2127    }
2128  }
2129
2130  /*!
2131    Send area from client to connected client(s)
2132  */
2133  void CDomain::sendArea()
2134  {
2135    if (!hasArea) return;
2136
2137    int ns, n, i, j, ind, nv, idx;
2138    std::list<CContextClient*>::iterator it;
2139
2140    for (it=clients.begin(); it!=clients.end(); ++it)
2141    {
2142      CContextClient* client = *it;
2143      int serverSize = client->serverSize;
2144
2145      // send area for each connected server
2146      CEventClient eventArea(getType(), EVENT_ID_AREA);
2147
2148      list<CMessage> list_msgsArea;
2149      list<CArray<double,1> > list_area;
2150
2151      std::unordered_map<int, vector<size_t> >::const_iterator it, iteMap;
2152      iteMap = indSrv_[serverSize].end();
2153      for (int k = 0; k < connectedServerRank_[serverSize].size(); ++k)
2154      {
2155        int nbData = 0;
2156        int rank = connectedServerRank_[serverSize][k];
2157        it = indSrv_[serverSize].find(rank);
2158        if (iteMap != it)
2159          nbData = it->second.size();
2160        list_area.push_back(CArray<double,1>(nbData));
2161
2162        const std::vector<size_t>& temp = it->second;
2163        for (n = 0; n < nbData; ++n)
2164        {
2165          idx = static_cast<int>(it->second[n]);
2166          list_area.back()(n) = areavalue(globalLocalIndexMap_[idx]);
2167        }
2168
2169        list_msgsArea.push_back(CMessage());
2170        list_msgsArea.back() << this->getId() << hasArea;
2171        list_msgsArea.back() << list_area.back();
2172        eventArea.push(rank, nbSenders[serverSize][rank], list_msgsArea.back());
2173      }
2174      client->sendEvent(eventArea);
2175    }
2176  }
2177
2178  /*!
2179    Send longitude and latitude from client to servers
2180    Each client send long and lat information to corresponding connected clients(s).
2181    Because longitude and latitude are optional, this function only called if latitude and longitude exist
2182  */
2183  void CDomain::sendLonLat()
2184  {
2185    if (!hasLonLat) return;
2186
2187    int ns, n, i, j, ind, nv, idx;
2188    std::list<CContextClient*>::iterator it;
2189    for (it=clients.begin(); it!=clients.end(); ++it)
2190    {
2191      CContextClient* client = *it;
2192      int serverSize = client->serverSize;
2193
2194      // send lon lat for each connected server
2195      CEventClient eventLon(getType(), EVENT_ID_LON);
2196      CEventClient eventLat(getType(), EVENT_ID_LAT);
2197
2198      list<CMessage> list_msgsLon, list_msgsLat;
2199      list<CArray<double,1> > list_lon, list_lat;
2200      list<CArray<double,2> > list_boundslon, list_boundslat;
2201
2202      std::unordered_map<int, vector<size_t> >::const_iterator it, iteMap;
2203      iteMap = indSrv_[serverSize].end();
2204      for (int k = 0; k < connectedServerRank_[serverSize].size(); ++k)
2205      {
2206        int nbData = 0;
2207        int rank = connectedServerRank_[serverSize][k];
2208        it = indSrv_[serverSize].find(rank);
2209        if (iteMap != it)
2210          nbData = it->second.size();
2211
2212        list_lon.push_back(CArray<double,1>(nbData));
2213        list_lat.push_back(CArray<double,1>(nbData));
2214
2215        if (hasBounds)
2216        {
2217          list_boundslon.push_back(CArray<double,2>(nvertex, nbData));
2218          list_boundslat.push_back(CArray<double,2>(nvertex, nbData));
2219        }
2220
2221        CArray<double,1>& lon = list_lon.back();
2222        CArray<double,1>& lat = list_lat.back();
2223        const std::vector<size_t>& temp = it->second;
2224        for (n = 0; n < nbData; ++n)
2225        {
2226          idx = static_cast<int>(it->second[n]);
2227          int localInd = globalLocalIndexMap_[idx];
2228          lon(n) = lonvalue(localInd);
2229          lat(n) = latvalue(localInd);
2230
2231          if (hasBounds)
2232          {
2233            CArray<double,2>& boundslon = list_boundslon.back();
2234            CArray<double,2>& boundslat = list_boundslat.back();
2235
2236            for (nv = 0; nv < nvertex; ++nv)
2237            {
2238              boundslon(nv, n) = bounds_lonvalue(nv, localInd);
2239              boundslat(nv, n) = bounds_latvalue(nv, localInd);
2240            }
2241          }
2242        }
2243
2244        list_msgsLon.push_back(CMessage());
2245        list_msgsLat.push_back(CMessage());
2246
2247        list_msgsLon.back() << this->getId() << hasLonLat;
2248        if (hasLonLat) 
2249          list_msgsLon.back() << list_lon.back();
2250        list_msgsLon.back()  << hasBounds;
2251        if (hasBounds)
2252        {
2253          list_msgsLon.back() << list_boundslon.back();
2254        }
2255
2256        list_msgsLat.back() << this->getId() << hasLonLat;
2257        if (hasLonLat)
2258          list_msgsLat.back() << list_lat.back();
2259        list_msgsLat.back() << hasBounds;
2260        if (hasBounds)
2261        {         
2262          list_msgsLat.back() << list_boundslat.back();
2263        }
2264
2265        eventLon.push(rank, nbSenders[serverSize][rank], list_msgsLon.back());
2266        eventLat.push(rank, nbSenders[serverSize][rank], list_msgsLat.back());
2267      }
2268      client->sendEvent(eventLon);
2269      client->sendEvent(eventLat);
2270    }
2271  }
2272
2273  /*!
2274    Send data index to corresponding connected clients.
2275    Data index can be compressed however, we always send decompressed data index
2276    and they will be compressed on receiving.
2277    The compressed index are represented with 1 and others are represented with -1
2278  */
2279  void CDomain::sendDataIndex()
2280  {
2281    int ns, n, i, j, ind, nv, idx;
2282    std::list<CContextClient*>::iterator it;
2283    for (it=clients.begin(); it!=clients.end(); ++it)
2284    {
2285      CContextClient* client = *it;
2286
2287      int serverSize = client->serverSize;
2288
2289      // send area for each connected server
2290      CEventClient eventDataIndex(getType(), EVENT_ID_DATA_INDEX);
2291
2292      list<CMessage> list_msgsDataIndex;
2293      list<CArray<int,1> > list_data_i_index, list_data_j_index;
2294
2295      int nbIndex = i_index.numElements();
2296      int niByIndex = max(i_index) - min(i_index) + 1;
2297      int njByIndex = max(j_index) - min(j_index) + 1; 
2298      int dataIindexBound = (1 == data_dim) ? (niByIndex * njByIndex) : niByIndex;
2299      int dataJindexBound = (1 == data_dim) ? (niByIndex * njByIndex) : njByIndex;
2300
2301     
2302      CArray<int,1> dataIIndex(nbIndex), dataJIndex(nbIndex);
2303      dataIIndex = -1; 
2304      dataJIndex = -1;
2305      ind = 0;
2306
2307      for (idx = 0; idx < data_i_index.numElements(); ++idx)
2308      {
2309        int dataIidx = data_i_index(idx) + data_ibegin;
2310        int dataJidx = data_j_index(idx) + data_jbegin;
2311        if ((0 <= dataIidx) && (dataIidx < dataIindexBound) &&
2312            (0 <= dataJidx) && (dataJidx < dataJindexBound))
2313        {
2314          dataIIndex((1 == data_dim) ? dataIidx : dataJidx * ni + dataIidx) = 1; //i_index(dataIidx);//dataIidx;
2315          dataJIndex((1 == data_dim) ? dataIidx : dataJidx * ni + dataIidx) = 1; //j_index(dataJidx);//         
2316        }
2317      }
2318
2319      std::unordered_map<int, vector<size_t> >::const_iterator it, iteMap;
2320      iteMap = indSrv_[serverSize].end();
2321      for (int k = 0; k < connectedServerRank_[serverSize].size(); ++k)
2322      {
2323        int nbData = 0;
2324        int rank = connectedServerRank_[serverSize][k];
2325        it = indSrv_[serverSize].find(rank);
2326        if (iteMap != it)
2327          nbData = it->second.size();
2328        list_data_i_index.push_back(CArray<int,1>(nbData));
2329        list_data_j_index.push_back(CArray<int,1>(nbData));
2330
2331        const std::vector<size_t>& temp = it->second;
2332        for (n = 0; n < nbData; ++n)
2333        {
2334          idx = static_cast<int>(it->second[n]);
2335          i = globalLocalIndexMap_[idx];
2336          list_data_i_index.back()(n) = dataIIndex(i);
2337          list_data_j_index.back()(n) = dataJIndex(i);
2338        }
2339
2340        list_msgsDataIndex.push_back(CMessage());
2341        list_msgsDataIndex.back() << this->getId();
2342        list_msgsDataIndex.back() << list_data_i_index.back() << list_data_j_index.back();
2343        eventDataIndex.push(rank, nbSenders[serverSize][rank], list_msgsDataIndex.back());
2344      }
2345      client->sendEvent(eventDataIndex);
2346    }
2347  }
2348 
2349  bool CDomain::dispatchEvent(CEventServer& event)
2350  {
2351    if (SuperClass::dispatchEvent(event)) return true;
2352    else
2353    {
2354      switch(event.type)
2355      {
2356        case EVENT_ID_SERVER_ATTRIBUT:
2357          recvDistributionAttributes(event);
2358          return true;
2359          break;
2360        case EVENT_ID_INDEX:
2361          recvIndex(event);
2362          return true;
2363          break;
2364        case EVENT_ID_MASK:
2365          recvMask(event);
2366          return true;
2367          break;
2368        case EVENT_ID_LON:
2369          recvLon(event);
2370          return true;
2371          break;
2372        case EVENT_ID_LAT:
2373          recvLat(event);
2374          return true;
2375          break;
2376        case EVENT_ID_AREA:
2377          recvArea(event);
2378          return true;
2379          break; 
2380        case EVENT_ID_DATA_INDEX:
2381          recvDataIndex(event);
2382          return true;
2383          break;
2384        default:
2385          ERROR("bool CDomain::dispatchEvent(CEventServer& event)",
2386                << "Unknown Event");
2387          return false;
2388       }
2389    }
2390  }
2391
2392  /*!
2393    Receive index event from clients(s)
2394    \param[in] event event contain info about rank and associated index
2395  */
2396  void CDomain::recvIndex(CEventServer& event)
2397  {
2398    string domainId;
2399    std::map<int, CBufferIn*> rankBuffers;
2400
2401    list<CEventServer::SSubEvent>::iterator it;
2402    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2403    {     
2404      CBufferIn* buffer = it->buffer;
2405      *buffer >> domainId;
2406      rankBuffers[it->rank] = buffer;       
2407    }
2408    get(domainId)->recvIndex(rankBuffers);
2409  }
2410
2411  /*!
2412    Receive index information from client(s). We use the global index for mapping index between
2413    sending clients and receiving clients.
2414    \param[in] rankBuffers rank of sending client and the corresponding receive buffer 
2415  */
2416  void CDomain::recvIndex(std::map<int, CBufferIn*>& rankBuffers)
2417  {
2418    int nbReceived = rankBuffers.size(), i, ind, index, type_int, iIndex, jIndex;
2419    recvClientRanks_.resize(nbReceived);       
2420
2421    std::map<int, CBufferIn*>::iterator it = rankBuffers.begin(), ite = rankBuffers.end();
2422    ind = 0;
2423    for (ind = 0; it != ite; ++it, ++ind)
2424    {       
2425       recvClientRanks_[ind] = it->first;
2426       CBufferIn& buffer = *(it->second);
2427       buffer >> type_int >> isCurvilinear >> indGlob_[it->first]; 
2428       type.setValue((type_attr::t_enum)type_int); // probleme des type enum avec les buffers : ToFix
2429    }
2430    int nbIndGlob = 0;
2431    for (i = 0; i < nbReceived; ++i)
2432    {
2433      nbIndGlob += indGlob_[recvClientRanks_[i]].numElements();
2434    }
2435   
2436    globalLocalIndexMap_.rehash(std::ceil(nbIndGlob/globalLocalIndexMap_.max_load_factor()));
2437    i_index.resize(nbIndGlob);
2438    j_index.resize(nbIndGlob);   
2439    int nbIndexGlobMax = nbIndGlob, nbIndLoc;
2440
2441    nbIndGlob = 0;
2442    for (i = 0; i < nbReceived; ++i)
2443    {
2444      CArray<int,1>& tmp = indGlob_[recvClientRanks_[i]];
2445      for (ind = 0; ind < tmp.numElements(); ++ind)
2446      {
2447         index = tmp(ind);
2448         if (0 == globalLocalIndexMap_.count(index))
2449         {
2450           iIndex = (index%ni_glo)-ibegin;
2451           iIndex = (iIndex < 0) ? 0 : iIndex;
2452           jIndex = (index/ni_glo)-jbegin;
2453           jIndex = (jIndex < 0) ? 0 : jIndex;
2454           nbIndLoc = iIndex + ni * jIndex;
2455           if (nbIndLoc < nbIndexGlobMax)
2456           {
2457             i_index(nbIndLoc) = index % ni_glo;
2458             j_index(nbIndLoc) = index / ni_glo;
2459             globalLocalIndexMap_[index] = nbIndLoc; 
2460             ++nbIndGlob;
2461           }
2462           // i_index(nbIndGlob) = index % ni_glo;
2463           // j_index(nbIndGlob) = index / ni_glo;
2464           // globalLocalIndexMap_[index] = nbIndGlob; 
2465           // ++nbIndGlob;
2466         } 
2467      } 
2468    } 
2469
2470    if (nbIndGlob==0)
2471    {
2472      i_index.resize(nbIndGlob);
2473      j_index.resize(nbIndGlob);
2474    }
2475    else
2476    {
2477      i_index.resizeAndPreserve(nbIndGlob);
2478      j_index.resizeAndPreserve(nbIndGlob);
2479    }
2480  }
2481
2482  /*!
2483    Receive attributes event from clients(s)
2484    \param[in] event event contain info about rank and associated attributes
2485  */
2486  void CDomain::recvDistributionAttributes(CEventServer& event)
2487  {
2488    CBufferIn* buffer=event.subEvents.begin()->buffer;
2489    string domainId ;
2490    *buffer>>domainId ;
2491    get(domainId)->recvDistributionAttributes(*buffer);
2492  }
2493
2494  /*!
2495    Receive attributes from client(s)
2496    \param[in] rank rank of client source
2497    \param[in] buffer message containing attributes info
2498  */
2499  void CDomain::recvDistributionAttributes(CBufferIn& buffer)
2500  {
2501    int ni_tmp, ibegin_tmp, nj_tmp, jbegin_tmp;
2502    int ni_glo_tmp, nj_glo_tmp;
2503    buffer >> isUnstructed_ >> ni_tmp >> ibegin_tmp >> nj_tmp >> jbegin_tmp
2504           >> ni_glo_tmp >> nj_glo_tmp
2505           >> isCompressible_;
2506
2507    ni.setValue(ni_tmp);
2508    ibegin.setValue(ibegin_tmp);
2509    nj.setValue(nj_tmp);
2510    jbegin.setValue(jbegin_tmp);
2511    ni_glo.setValue(ni_glo_tmp);
2512    nj_glo.setValue(nj_glo_tmp);
2513
2514  }
2515
2516  /*!
2517    Receive area event from clients(s)
2518    \param[in] event event contain info about rank and associated area
2519  */
2520  void CDomain::recvMask(CEventServer& event)
2521  {
2522    string domainId;
2523    std::map<int, CBufferIn*> rankBuffers;
2524
2525    list<CEventServer::SSubEvent>::iterator it;
2526    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2527    {     
2528      CBufferIn* buffer = it->buffer;
2529      *buffer >> domainId;
2530      rankBuffers[it->rank] = buffer;     
2531    }
2532    get(domainId)->recvMask(rankBuffers);
2533  }
2534
2535
2536  /*!
2537    Receive mask information from client(s)
2538    \param[in] rankBuffers rank of sending client and the corresponding receive buffer 
2539  */
2540  void CDomain::recvMask(std::map<int, CBufferIn*>& rankBuffers)
2541  {
2542    int nbReceived = rankBuffers.size(), i, ind, index, lInd;
2543    if (nbReceived != recvClientRanks_.size())
2544      ERROR("void CDomain::recvMask(std::map<int, CBufferIn*>& rankBuffers)",
2545           << "The number of sending clients is not correct."
2546           << "Expected number: " << recvClientRanks_.size() << " but received " << nbReceived);
2547
2548    vector<CArray<bool,1> > recvMaskValue(nbReceived);     
2549    for (i = 0; i < recvClientRanks_.size(); ++i)
2550    {
2551      int rank = recvClientRanks_[i];
2552      CBufferIn& buffer = *(rankBuffers[rank]);     
2553      buffer >> recvMaskValue[i];
2554    }
2555
2556    int nbMaskInd = 0;
2557    for (i = 0; i < nbReceived; ++i)
2558    {
2559      nbMaskInd += recvMaskValue[i].numElements();
2560    }
2561 
2562    if (nbMaskInd != globalLocalIndexMap_.size())
2563      info (0) << "If domain " << this->getDomainOutputName() <<" does not have overlapped regions between processes "
2564               << "something must be wrong with mask index "<< std::endl;
2565
2566    nbMaskInd = globalLocalIndexMap_.size();
2567    mask_1d.resize(nbMaskInd);
2568    domainMask.resize(nbMaskInd);
2569    mask_1d = false;
2570   
2571    for (i = 0; i < nbReceived; ++i)
2572    {
2573      CArray<int,1>& tmpInd = indGlob_[recvClientRanks_[i]];
2574      CArray<bool,1>& tmp = recvMaskValue[i];
2575      for (ind = 0; ind < tmp.numElements(); ++ind)
2576      {
2577        lInd = globalLocalIndexMap_[size_t(tmpInd(ind))];
2578        if (!mask_1d(lInd)) // Only rewrite mask_1d if it's not true
2579          mask_1d(lInd) = tmp(ind);
2580      }
2581    }
2582    domainMask=mask_1d ;
2583  }
2584
2585  /*!
2586    Receive longitude event from clients(s)
2587    \param[in] event event contain info about rank and associated longitude
2588  */
2589  void CDomain::recvLon(CEventServer& event)
2590  {
2591    string domainId;
2592    std::map<int, CBufferIn*> rankBuffers;
2593
2594    list<CEventServer::SSubEvent>::iterator it;
2595    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2596    {     
2597      CBufferIn* buffer = it->buffer;
2598      *buffer >> domainId;
2599      rankBuffers[it->rank] = buffer;       
2600    }
2601    get(domainId)->recvLon(rankBuffers);
2602  }
2603
2604  /*!
2605    Receive longitude information from client(s)
2606    \param[in] rankBuffers rank of sending client and the corresponding receive buffer 
2607  */
2608  void CDomain::recvLon(std::map<int, CBufferIn*>& rankBuffers)
2609  {
2610    int nbReceived = rankBuffers.size(), i, ind, index, iindex, jindex, lInd;
2611    if (nbReceived != recvClientRanks_.size())
2612      ERROR("void CDomain::recvLon(std::map<int, CBufferIn*>& rankBuffers)",
2613           << "The number of sending clients is not correct."
2614           << "Expected number: " << recvClientRanks_.size() << " but received " << nbReceived);
2615
2616    vector<CArray<double,1> > recvLonValue(nbReceived);
2617    vector<CArray<double,2> > recvBoundsLonValue(nbReceived);   
2618    for (i = 0; i < recvClientRanks_.size(); ++i)
2619    {
2620      int rank = recvClientRanks_[i];
2621      CBufferIn& buffer = *(rankBuffers[rank]);
2622      buffer >> hasLonLat;
2623      if (hasLonLat)
2624        buffer >> recvLonValue[i];
2625      buffer >> hasBounds;
2626      if (hasBounds)
2627        buffer >> recvBoundsLonValue[i];
2628    }
2629
2630    if (hasLonLat)
2631    {
2632      int nbLonInd = 0;
2633      for (i = 0; i < nbReceived; ++i)
2634      {
2635        nbLonInd += recvLonValue[i].numElements();
2636      }
2637   
2638      if (nbLonInd != globalLocalIndexMap_.size())
2639        info (0) << "If domain " << this->getDomainOutputName() <<" does not have overlapped regions between processes "
2640                 << "something must be wrong with longitude index "<< std::endl;
2641
2642      nbLonInd = globalLocalIndexMap_.size();
2643      lonvalue.resize(nbLonInd);
2644      if (hasBounds)
2645      {
2646        bounds_lonvalue.resize(nvertex,nbLonInd);
2647        bounds_lonvalue = 0.;
2648      }
2649
2650      nbLonInd = 0;
2651      for (i = 0; i < nbReceived; ++i)
2652      {
2653        CArray<int,1>& tmpInd = indGlob_[recvClientRanks_[i]];
2654        CArray<double,1>& tmp = recvLonValue[i];
2655        for (ind = 0; ind < tmp.numElements(); ++ind)
2656        {
2657          lInd = globalLocalIndexMap_[size_t(tmpInd(ind))];
2658          lonvalue(lInd) = tmp(ind); 
2659           if (hasBounds)
2660           {         
2661            for (int nv = 0; nv < nvertex; ++nv)
2662              bounds_lonvalue(nv, lInd) = recvBoundsLonValue[i](nv, ind);
2663           }                 
2664        }
2665      }       
2666    }
2667  }
2668
2669  /*!
2670    Receive latitude event from clients(s)
2671    \param[in] event event contain info about rank and associated latitude
2672  */
2673  void CDomain::recvLat(CEventServer& event)
2674  {
2675    string domainId;
2676    std::map<int, CBufferIn*> rankBuffers;
2677
2678    list<CEventServer::SSubEvent>::iterator it;
2679    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2680    {     
2681      CBufferIn* buffer = it->buffer;
2682      *buffer >> domainId;
2683      rankBuffers[it->rank] = buffer;   
2684    }
2685    get(domainId)->recvLat(rankBuffers);
2686  }
2687
2688  /*!
2689    Receive latitude information from client(s)
2690    \param[in] rankBuffers rank of sending client and the corresponding receive buffer 
2691  */
2692  void CDomain::recvLat(std::map<int, CBufferIn*>& rankBuffers)
2693  {
2694    int nbReceived = rankBuffers.size(), i, ind, index, iindex, jindex, lInd;
2695    if (nbReceived != recvClientRanks_.size())
2696      ERROR("void CDomain::recvLat(std::map<int, CBufferIn*>& rankBuffers)",
2697           << "The number of sending clients is not correct."
2698           << "Expected number: " << recvClientRanks_.size() << " but received " << nbReceived);
2699
2700    vector<CArray<double,1> > recvLatValue(nbReceived);
2701    vector<CArray<double,2> > recvBoundsLatValue(nbReceived);   
2702    for (i = 0; i < recvClientRanks_.size(); ++i)
2703    {
2704      int rank = recvClientRanks_[i];
2705      CBufferIn& buffer = *(rankBuffers[rank]);
2706      buffer >> hasLonLat;
2707      if (hasLonLat)
2708        buffer >> recvLatValue[i];
2709      buffer >> hasBounds;
2710      if (hasBounds)
2711        buffer >> recvBoundsLatValue[i];
2712    }
2713
2714    if (hasLonLat)
2715    {
2716      int nbLatInd = 0;
2717      for (i = 0; i < nbReceived; ++i)
2718      {
2719        nbLatInd += recvLatValue[i].numElements();
2720      }
2721   
2722      if (nbLatInd != globalLocalIndexMap_.size())
2723        info (0) << "If domain " << this->getDomainOutputName() <<" does not have overlapped regions between processes "
2724                << "something must be wrong with latitude index "<< std::endl;
2725
2726      nbLatInd = globalLocalIndexMap_.size();
2727      latvalue.resize(nbLatInd);
2728      if (hasBounds)
2729      {
2730        bounds_latvalue.resize(nvertex,nbLatInd);
2731        bounds_latvalue = 0. ;
2732      }
2733
2734      nbLatInd = 0;
2735      for (i = 0; i < nbReceived; ++i)
2736      {
2737        CArray<int,1>& tmpInd = indGlob_[recvClientRanks_[i]];
2738        CArray<double,1>& tmp = recvLatValue[i];
2739        for (ind = 0; ind < tmp.numElements(); ++ind)
2740        {
2741          lInd = globalLocalIndexMap_[size_t(tmpInd(ind))];
2742          latvalue(lInd) = tmp(ind);   
2743           if (hasBounds)
2744           {
2745            CArray<double,2>& boundslat = recvBoundsLatValue[i];
2746            for (int nv = 0; nv < nvertex; ++nv)
2747              bounds_latvalue(nv, lInd) = boundslat(nv, ind);
2748           }   
2749          ++nbLatInd;
2750        }
2751      }       
2752    }
2753  }
2754
2755  /*!
2756    Receive area event from clients(s)
2757    \param[in] event event contain info about rank and associated area
2758  */
2759  void CDomain::recvArea(CEventServer& event)
2760  {
2761    string domainId;
2762    std::map<int, CBufferIn*> rankBuffers;
2763
2764    list<CEventServer::SSubEvent>::iterator it;
2765    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2766    {     
2767      CBufferIn* buffer = it->buffer;
2768      *buffer >> domainId;
2769      rankBuffers[it->rank] = buffer;     
2770    }
2771    get(domainId)->recvArea(rankBuffers);
2772  }
2773
2774  /*!
2775    Receive area information from client(s)
2776    \param[in] rankBuffers rank of sending client and the corresponding receive buffer     
2777  */
2778  void CDomain::recvArea(std::map<int, CBufferIn*>& rankBuffers)
2779  {
2780    int nbReceived = rankBuffers.size(), i, ind, index, lInd;
2781    if (nbReceived != recvClientRanks_.size())
2782      ERROR("void CDomain::recvArea(std::map<int, CBufferIn*>& rankBuffers)",
2783           << "The number of sending clients is not correct."
2784           << "Expected number: " << recvClientRanks_.size() << " but received " << nbReceived);
2785
2786    vector<CArray<double,1> > recvAreaValue(nbReceived);     
2787    for (i = 0; i < recvClientRanks_.size(); ++i)
2788    {
2789      int rank = recvClientRanks_[i];
2790      CBufferIn& buffer = *(rankBuffers[rank]);     
2791      buffer >> hasArea;
2792      if (hasArea)
2793        buffer >> recvAreaValue[i];
2794    }
2795
2796    if (hasArea)
2797    {
2798      int nbAreaInd = 0;
2799      for (i = 0; i < nbReceived; ++i)
2800      {     
2801        nbAreaInd += recvAreaValue[i].numElements();
2802      }
2803
2804      if (nbAreaInd != globalLocalIndexMap_.size())
2805        info (0) << "If domain " << this->getDomainOutputName() <<" does not have overlapped regions between processes "
2806                 << "something must be wrong with area index "<< std::endl;
2807
2808      nbAreaInd = globalLocalIndexMap_.size();
2809      areavalue.resize(nbAreaInd);
2810      nbAreaInd = 0;     
2811      for (i = 0; i < nbReceived; ++i)
2812      {
2813        CArray<int,1>& tmpInd = indGlob_[recvClientRanks_[i]];
2814        CArray<double,1>& tmp = recvAreaValue[i];
2815        for (ind = 0; ind < tmp.numElements(); ++ind)
2816        {
2817          lInd = globalLocalIndexMap_[size_t(tmpInd(ind))];
2818          areavalue(lInd) = tmp(ind);         
2819        }
2820      }
2821     
2822    }
2823  }
2824
2825  /*!
2826    Compare two domain objects.
2827    They are equal if only if they have identical attributes as well as their values.
2828    Moreover, they must have the same transformations.
2829  \param [in] domain Compared domain
2830  \return result of the comparison
2831  */
2832  bool CDomain::isEqual(CDomain* obj)
2833  {
2834    vector<StdString> excludedAttr;
2835    excludedAttr.push_back("domain_ref");
2836    bool objEqual = SuperClass::isEqual(obj, excludedAttr);
2837    if (!objEqual) return objEqual;
2838
2839    TransMapTypes thisTrans = this->getAllTransformations();
2840    TransMapTypes objTrans  = obj->getAllTransformations();
2841
2842    TransMapTypes::const_iterator it, itb, ite;
2843    std::vector<ETranformationType> thisTransType, objTransType;
2844    for (it = thisTrans.begin(); it != thisTrans.end(); ++it)
2845      thisTransType.push_back(it->first);
2846    for (it = objTrans.begin(); it != objTrans.end(); ++it)
2847      objTransType.push_back(it->first);
2848
2849    if (thisTransType.size() != objTransType.size()) return false;
2850    for (int idx = 0; idx < thisTransType.size(); ++idx)
2851      objEqual &= (thisTransType[idx] == objTransType[idx]);
2852
2853    return objEqual;
2854  }
2855
2856  /*!
2857    Receive data index event from clients(s)
2858    \param[in] event event contain info about rank and associated index
2859  */
2860  void CDomain::recvDataIndex(CEventServer& event)
2861  {
2862    string domainId;
2863    std::map<int, CBufferIn*> rankBuffers;
2864
2865    list<CEventServer::SSubEvent>::iterator it;
2866    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
2867    {     
2868      CBufferIn* buffer = it->buffer;
2869      *buffer >> domainId;
2870      rankBuffers[it->rank] = buffer;       
2871    }
2872    get(domainId)->recvDataIndex(rankBuffers);
2873  }
2874
2875  /*!
2876    Receive data index information from client(s)
2877    A client receives data index from different clients to rebuild its own data index.
2878    Because we use global index + mask info to calculate the sending data to client(s),
2879    this data index must be updated with mask info (maybe it will change in the future)
2880    Because the data index is local, to rebuild data index of received client, we should use global index along with.
2881
2882    \param[in] rankBuffers rank of sending client and the corresponding receive buffer     
2883  */
2884  void CDomain::recvDataIndex(std::map<int, CBufferIn*>& rankBuffers)
2885  {
2886    int nbReceived = rankBuffers.size(), i, ind, index, indexI, indexJ, type_int, lInd;   
2887    if (nbReceived != recvClientRanks_.size())
2888      ERROR("void CDomain::recvDataIndex(std::map<int, CBufferIn*>& rankBuffers)",
2889           << "The number of sending clients is not correct."
2890           << "Expected number: " << recvClientRanks_.size() << " but received " << nbReceived);
2891
2892    vector<CArray<int,1> > recvDataIIndex(nbReceived),recvDataJIndex(nbReceived);     
2893    for (i = 0; i < recvClientRanks_.size(); ++i)
2894    {
2895      int rank = recvClientRanks_[i];
2896      CBufferIn& buffer = *(rankBuffers[rank]);
2897      buffer >> recvDataIIndex[i];
2898      buffer >> recvDataJIndex[i];
2899    }
2900   
2901    int nbIndex = i_index.numElements();
2902    CArray<int,1> dataIIndex(nbIndex), dataJIndex(nbIndex);
2903    dataIIndex = -1; dataJIndex = -1;
2904     
2905    nbIndex = 0;
2906    for (i = 0; i < nbReceived; ++i)
2907    {     
2908      CArray<int,1>& tmpInd = indGlob_[recvClientRanks_[i]];
2909      CArray<int,1>& tmpI = recvDataIIndex[i];   
2910      CArray<int,1>& tmpJ = recvDataJIndex[i];     
2911      if ((tmpI.numElements() != tmpInd.numElements()) || (tmpJ.numElements() != tmpInd.numElements()))
2912          ERROR("void CDomain::recvDataIndex(std::map<int, CBufferIn*>& rankBuffers)",
2913             << "The number of global received index is not coherent with the number of received data index."
2914             << "Expected number of global index: " << tmpI.numElements() << " but received " << tmpInd.numElements());
2915
2916      for (ind = 0; ind < tmpI.numElements(); ++ind)
2917      {
2918         lInd = globalLocalIndexMap_[size_t(tmpInd(ind))];
2919         dataIIndex(lInd) = (-1 == dataIIndex(lInd)) ? tmpI(ind) : dataIIndex(lInd); // Only fill in dataIndex if there is no data
2920         dataJIndex(lInd) = (-1 == dataJIndex(lInd)) ? tmpJ(ind) : dataJIndex(lInd); 
2921
2922         if (!domainMask(lInd))   // Include mask info into data index on the RECEIVE getServerDimensionSizes
2923         {
2924           dataIIndex(lInd) = dataJIndex(lInd) = -1;
2925         }
2926      } 
2927    }
2928
2929    int nbCompressedData = 0; 
2930    for (ind = 0; ind < dataIIndex.numElements(); ++ind)
2931    {
2932       indexI = dataIIndex(ind); indexJ = dataJIndex(ind);
2933       if ((0 <= indexI) && (0 <= indexJ))
2934         ++nbCompressedData;
2935    }       
2936 
2937    data_i_index.resize(nbCompressedData);
2938    data_j_index.resize(nbCompressedData);
2939
2940    nbCompressedData = 0; 
2941    for (ind = 0; ind < dataIIndex.numElements(); ++ind)
2942    {
2943       indexI = dataIIndex(ind); indexJ = dataJIndex(ind);
2944       if ((0 <= indexI) && (0 <= indexJ))
2945       {
2946          data_i_index(nbCompressedData) = (1 == data_dim) ? ind : ind % ni;
2947          data_j_index(nbCompressedData) = (1 == data_dim) ? 0   : ind / ni;
2948         ++nbCompressedData;
2949       }
2950    }
2951
2952    // Reset data_ibegin, data_jbegin
2953    data_ibegin.setValue(0);
2954    data_jbegin.setValue(0);
2955  }
2956
2957  CTransformation<CDomain>* CDomain::addTransformation(ETranformationType transType, const StdString& id)
2958  {
2959    transformationMap_.push_back(std::make_pair(transType, CTransformation<CDomain>::createTransformation(transType,id)));
2960    return transformationMap_.back().second;
2961  }
2962
2963  /*!
2964    Check whether a domain has transformation
2965    \return true if domain has transformation
2966  */
2967  bool CDomain::hasTransformation()
2968  {
2969    return (!transformationMap_.empty());
2970  }
2971
2972  /*!
2973    Set transformation for current domain. It's the method to move transformation in hierarchy
2974    \param [in] domTrans transformation on domain
2975  */
2976  void CDomain::setTransformations(const TransMapTypes& domTrans)
2977  {
2978    transformationMap_ = domTrans;
2979  }
2980
2981  /*!
2982    Get all transformation current domain has
2983    \return all transformation
2984  */
2985  CDomain::TransMapTypes CDomain::getAllTransformations(void)
2986  {
2987    return transformationMap_;
2988  }
2989
2990  void CDomain::duplicateTransformation(CDomain* src)
2991  {
2992    if (src->hasTransformation())
2993    {
2994      this->setTransformations(src->getAllTransformations());
2995    }
2996  }
2997
2998  /*!
2999   * Go through the hierarchy to find the domain from which the transformations must be inherited
3000   */
3001  void CDomain::solveInheritanceTransformation()
3002  {
3003    if (hasTransformation() || !hasDirectDomainReference())
3004      return;
3005
3006    CDomain* domain = this;
3007    std::vector<CDomain*> refDomains;
3008    while (!domain->hasTransformation() && domain->hasDirectDomainReference())
3009    {
3010      refDomains.push_back(domain);
3011      domain = domain->getDirectDomainReference();
3012    }
3013
3014    if (domain->hasTransformation())
3015      for (size_t i = 0; i < refDomains.size(); ++i)
3016        refDomains[i]->setTransformations(domain->getAllTransformations());
3017  }
3018
3019  void CDomain::setContextClient(CContextClient* contextClient)
3020  {
3021    if (clientsSet.find(contextClient)==clientsSet.end())
3022    {
3023      clients.push_back(contextClient) ;
3024      clientsSet.insert(contextClient);
3025    }
3026  }
3027
3028  /*!
3029    Parse children nodes of a domain in xml file.
3030    Whenver there is a new transformation, its type and name should be added into this function
3031    \param node child node to process
3032  */
3033  void CDomain::parse(xml::CXMLNode & node)
3034  {
3035    SuperClass::parse(node);
3036
3037    if (node.goToChildElement())
3038    {
3039      StdString nodeElementName;
3040      do
3041      {
3042        StdString nodeId("");
3043        if (node.getAttributes().end() != node.getAttributes().find("id"))
3044        { nodeId = node.getAttributes()["id"]; }
3045
3046        nodeElementName = node.getElementName();
3047        std::map<StdString, ETranformationType>::const_iterator ite = transformationMapList_.end(), it;
3048        it = transformationMapList_.find(nodeElementName);
3049        if (ite != it)
3050        {
3051          transformationMap_.push_back(std::make_pair(it->second, CTransformation<CDomain>::createTransformation(it->second,
3052                                                                                                                nodeId,
3053                                                                                                                &node)));
3054        }
3055        else
3056        {
3057          ERROR("void CDomain::parse(xml::CXMLNode & node)",
3058                << "The transformation " << nodeElementName << " has not been supported yet.");
3059        }
3060      } while (node.goToNextElement()) ;
3061      node.goToParentElement();
3062    }
3063  }
3064   //----------------------------------------------------------------
3065
3066   DEFINE_REF_FUNC(Domain,domain)
3067
3068   ///---------------------------------------------------------------
3069
3070} // namespace xios
Note: See TracBrowser for help on using the repository browser.