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

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

MARK: Dynamic workflow graph developement. Branch up to date with trunk @1676. Generic testcase not yet tested

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