source: XIOS/dev/dev_trunk_omp/src/node/domain.cpp @ 1619

Last change on this file since 1619 was 1601, checked in by yushan, 6 years ago

branch_openmp merged with trunk r1597

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