source: XIOS/dev/dev_ym/XIOS_COUPLING/src/node/context.cpp @ 1870

Last change on this file since 1870 was 1870, checked in by ymipsl, 4 years ago

Some update on XIOS_COUPLING branch...

YM

  • Property copyright set to
    Software name : XIOS (Xml I/O Server)
    http://forge.ipsl.jussieu.fr/ioserver
    Creation date : January 2009
    Licence : CeCCIL version2
    see license file in root directory : Licence_CeCILL_V2-en.txt
    or http://www.cecill.info/licences/Licence_CeCILL_V2-en.html
    Holder : CEA/LSCE (Laboratoire des Sciences du CLimat et de l'Environnement)
    CNRS/IPSL (Institut Pierre Simon Laplace)
    Project Manager : Yann Meurdesoif
    yann.meurdesoif@cea.fr
File size: 96.9 KB
Line 
1#include "context.hpp"
2#include "attribute_template.hpp"
3#include "object_template.hpp"
4#include "group_template.hpp"
5
6#include "calendar_type.hpp"
7#include "duration.hpp"
8
9#include "context_client.hpp"
10#include "context_server.hpp"
11#include "nc4_data_output.hpp"
12#include "node_type.hpp"
13#include "message.hpp"
14#include "type.hpp"
15#include "xios_spl.hpp"
16#include "timer.hpp"
17#include "memtrack.hpp"
18#include <limits>
19#include <fstream>
20#include "server.hpp"
21#include "distribute_file_server2.hpp"
22#include "services_manager.hpp"
23#include "contexts_manager.hpp"
24#include "cxios.hpp"
25#include "client.hpp"
26#include "coupler_in.hpp"
27#include "coupler_out.hpp"
28
29namespace xios {
30
31  std::shared_ptr<CContextGroup> CContext::root;
32
33   /// ////////////////////// Définitions ////////////////////// ///
34
35   CContext::CContext(void)
36      : CObjectTemplate<CContext>(), CContextAttributes()
37      , calendar(), hasClient(false), hasServer(false)
38      , isPostProcessed(false), finalized(false)
39      , client(nullptr), server(nullptr)
40      , allProcessed(false), countChildContextFinalized_(0), isProcessingEvent_(false)
41
42   { /* Ne rien faire de plus */ }
43
44   CContext::CContext(const StdString & id)
45      : CObjectTemplate<CContext>(id), CContextAttributes()
46      , calendar(), hasClient(false), hasServer(false)
47      , isPostProcessed(false), finalized(false)
48      , client(nullptr), server(nullptr)
49      , allProcessed(false), countChildContextFinalized_(0), isProcessingEvent_(false)
50   { /* Ne rien faire de plus */ }
51
52   CContext::~CContext(void)
53   {
54     delete client;
55     delete server;
56     for (std::vector<CContextClient*>::iterator it = clientPrimServer.begin(); it != clientPrimServer.end(); it++)  delete *it;
57     for (std::vector<CContextServer*>::iterator it = serverPrimServer.begin(); it != serverPrimServer.end(); it++)  delete *it;
58
59   }
60
61   //----------------------------------------------------------------
62   //! Get name of context
63   StdString CContext::GetName(void)   { return (StdString("context")); }
64   StdString CContext::GetDefName(void){ return (CContext::GetName()); }
65   ENodeType CContext::GetType(void)   { return (eContext); }
66
67   //----------------------------------------------------------------
68
69   /*!
70   \brief Get context group (context root)
71   \return Context root
72   */
73   CContextGroup* CContext::getRoot(void)
74   TRY
75   {
76      if (root.get()==NULL) root=std::shared_ptr<CContextGroup>(new CContextGroup(xml::CXMLNode::GetRootName()));
77      return root.get();
78   }
79   CATCH
80
81   //----------------------------------------------------------------
82
83   /*!
84   \brief Get calendar of a context
85   \return Calendar
86   */
87   std::shared_ptr<CCalendar> CContext::getCalendar(void) const
88   TRY
89   {
90      return (this->calendar);
91   }
92   CATCH
93
94   //----------------------------------------------------------------
95
96   /*!
97   \brief Set a context with a calendar
98   \param[in] newCalendar new calendar
99   */
100   void CContext::setCalendar(std::shared_ptr<CCalendar> newCalendar)
101   TRY
102   {
103      this->calendar = newCalendar;
104   }
105   CATCH_DUMP_ATTR
106
107   //----------------------------------------------------------------
108   /*!
109   \brief Parse xml file and write information into context object
110   \param [in] node xmld node corresponding in xml file
111   */
112   void CContext::parse(xml::CXMLNode & node)
113   TRY
114   {
115      CContext::SuperClass::parse(node);
116
117      // PARSING POUR GESTION DES ENFANTS
118      xml::THashAttributes attributes = node.getAttributes();
119
120      if (attributes.end() != attributes.find("src"))
121      {
122         StdIFStream ifs ( attributes["src"].c_str() , StdIFStream::in );
123         if ( (ifs.rdstate() & std::ifstream::failbit ) != 0 )
124            ERROR("void CContext::parse(xml::CXMLNode & node)",
125                  <<endl<< "Can not open <"<<attributes["src"].c_str()<<"> file" );
126         if (!ifs.good())
127            ERROR("CContext::parse(xml::CXMLNode & node)",
128                  << "[ filename = " << attributes["src"] << " ] Bad xml stream !");
129         xml::CXMLParser::ParseInclude(ifs, attributes["src"], *this);
130      }
131
132      if (node.getElementName().compare(CContext::GetName()))
133         DEBUG("Le noeud is wrong defined but will be considered as a context !");
134
135      if (!(node.goToChildElement()))
136      {
137         DEBUG("Le context ne contient pas d'enfant !");
138      }
139      else
140      {
141         do { // Parcours des contextes pour traitement.
142
143            StdString name = node.getElementName();
144            attributes.clear();
145            attributes = node.getAttributes();
146
147            if (attributes.end() != attributes.find("id"))
148            {
149              DEBUG(<< "Definition node has an id,"
150                    << "it will not be taking account !");
151            }
152
153#define DECLARE_NODE(Name_, name_)    \
154   if (name.compare(C##Name_##Definition::GetDefName()) == 0) \
155   { C##Name_##Definition::create(C##Name_##Definition::GetDefName()) -> parse(node); continue; }
156#define DECLARE_NODE_PAR(Name_, name_)
157#include "node_type.conf"
158
159            DEBUG(<< "The element \'"     << name
160                  << "\' in the context \'" << CContext::getCurrent()->getId()
161                  << "\' is not a definition !");
162
163         } while (node.goToNextElement());
164
165         node.goToParentElement(); // Retour au parent
166      }
167   }
168   CATCH_DUMP_ATTR
169
170   //----------------------------------------------------------------
171   //! Show tree structure of context
172   void CContext::ShowTree(StdOStream & out)
173   TRY
174   {
175      StdString currentContextId = CContext::getCurrent() -> getId();
176      std::vector<CContext*> def_vector =
177         CContext::getRoot()->getChildList();
178      std::vector<CContext*>::iterator
179         it = def_vector.begin(), end = def_vector.end();
180
181      out << "<? xml version=\"1.0\" ?>" << std::endl;
182      out << "<"  << xml::CXMLNode::GetRootName() << " >" << std::endl;
183
184      for (; it != end; it++)
185      {
186         CContext* context = *it;
187         CContext::setCurrent(context->getId());
188         out << *context << std::endl;
189      }
190
191      out << "</" << xml::CXMLNode::GetRootName() << " >" << std::endl;
192      CContext::setCurrent(currentContextId);
193   }
194   CATCH
195
196   //----------------------------------------------------------------
197
198   //! Convert context object into string (to print)
199   StdString CContext::toString(void) const
200   TRY
201   {
202      StdOStringStream oss;
203      oss << "<" << CContext::GetName()
204          << " id=\"" << this->getId() << "\" "
205          << SuperClassAttribute::toString() << ">" << std::endl;
206      if (!this->hasChild())
207      {
208         //oss << "<!-- No definition -->" << std::endl; // fait planter l'incrémentation
209      }
210      else
211      {
212
213#define DECLARE_NODE(Name_, name_)    \
214   if (C##Name_##Definition::has(C##Name_##Definition::GetDefName())) \
215   oss << * C##Name_##Definition::get(C##Name_##Definition::GetDefName()) << std::endl;
216#define DECLARE_NODE_PAR(Name_, name_)
217#include "node_type.conf"
218
219      }
220      oss << "</" << CContext::GetName() << " >";
221      return (oss.str());
222   }
223   CATCH
224
225   //----------------------------------------------------------------
226
227   /*!
228   \brief Find all inheritace among objects in a context.
229   \param [in] apply (true) write attributes of parent into ones of child if they are empty
230                     (false) write attributes of parent into a new container of child
231   \param [in] parent unused
232   */
233   void CContext::solveDescInheritance(bool apply, const CAttributeMap * const UNUSED(parent))
234   TRY
235   {
236#define DECLARE_NODE(Name_, name_)    \
237   if (C##Name_##Definition::has(C##Name_##Definition::GetDefName())) \
238     C##Name_##Definition::get(C##Name_##Definition::GetDefName())->solveDescInheritance(apply);
239#define DECLARE_NODE_PAR(Name_, name_)
240#include "node_type.conf"
241   }
242   CATCH_DUMP_ATTR
243
244   //----------------------------------------------------------------
245
246   //! Verify if all root definition in the context have child.
247   bool CContext::hasChild(void) const
248   TRY
249   {
250      return (
251#define DECLARE_NODE(Name_, name_)    \
252   C##Name_##Definition::has(C##Name_##Definition::GetDefName())   ||
253#define DECLARE_NODE_PAR(Name_, name_)
254#include "node_type.conf"
255      false);
256}
257   CATCH
258
259   //----------------------------------------------------------------
260
261   void CContext::CleanTree(void)
262   TRY
263   {
264#define DECLARE_NODE(Name_, name_) C##Name_##Definition::ClearAllAttributes();
265#define DECLARE_NODE_PAR(Name_, name_)
266#include "node_type.conf"
267   }
268   CATCH
269
270   ///---------------------------------------------------------------
271
272
273   void CContext::setClientServerBuffer(vector<CField*>& fields, bool bufferForWriting)
274   TRY
275   {
276      // Estimated minimum event size for small events (20 is an arbitrary constant just for safety)
277     const size_t minEventSize = CEventClient::headerSize + 20 * sizeof(int);
278      // Ensure there is at least some room for 20 of such events in the buffers
279     size_t minBufferSize = std::max(CXios::minBufferSize, 20 * minEventSize);
280
281#define DECLARE_NODE(Name_, name_)    \
282     if (minBufferSize < sizeof(C##Name_##Definition)) minBufferSize = sizeof(C##Name_##Definition);
283#define DECLARE_NODE_PAR(Name_, name_)
284#include "node_type.conf"
285#undef DECLARE_NODE
286#undef DECLARE_NODE_PAR
287
288
289     map<CContextClient*,map<int,size_t>> dataSize ;
290     map<CContextClient*,map<int,size_t>> maxEventSize ;
291     map<CContextClient*,map<int,size_t>> attributesSize ; 
292
293     for(auto field : fields)
294     {
295       field->setContextClientDataBufferSize(dataSize, maxEventSize, bufferForWriting) ;
296       field->setContextClientAttributesBufferSize(attributesSize, maxEventSize, bufferForWriting) ;
297     }
298     
299
300     for(auto& it : attributesSize)
301     {
302       auto contextClient = it.first ;
303       auto& contextDataSize =  dataSize[contextClient] ;
304       auto& contextAttributesSize =  attributesSize[contextClient] ;
305       auto& contextMaxEventSize =  maxEventSize[contextClient] ;
306   
307       for (auto& it : contextAttributesSize)
308       {
309         auto serverRank=it.first ;
310         auto& buffer = contextAttributesSize[serverRank] ;
311         if (contextDataSize[serverRank] > buffer) buffer=contextDataSize[serverRank] ;
312         buffer *= CXios::bufferSizeFactor;
313         if (buffer < minBufferSize) buffer = minBufferSize;
314         if (buffer > CXios::maxBufferSize ) buffer = CXios::maxBufferSize;
315       }
316
317       // Leaders will have to send some control events so ensure there is some room for those in the buffers
318       if (contextClient->isServerLeader())
319         for(auto& rank : contextClient->getRanksServerLeader())
320           if (!contextAttributesSize.count(rank))
321           {
322             contextAttributesSize[rank] = minBufferSize;
323             contextMaxEventSize[rank] = minEventSize;
324           }
325     
326       contextClient->setBufferSize(contextAttributesSize, contextMaxEventSize);   
327     }
328   }
329   CATCH_DUMP_ATTR
330
331
332    /*!
333    Sets client buffers.
334    \param [in] contextClient
335    \param [in] bufferForWriting True if buffers are used for sending data for writing
336    This flag is only true for client and server-1 for communication with server-2
337  */
338  // ym obsolete to be removed
339   void CContext::setClientServerBuffer(CContextClient* contextClient, bool bufferForWriting)
340   TRY
341   {
342      // Estimated minimum event size for small events (20 is an arbitrary constant just for safety)
343     const size_t minEventSize = CEventClient::headerSize + 20 * sizeof(int);
344
345      // Ensure there is at least some room for 20 of such events in the buffers
346      size_t minBufferSize = std::max(CXios::minBufferSize, 20 * minEventSize);
347
348#define DECLARE_NODE(Name_, name_)    \
349     if (minBufferSize < sizeof(C##Name_##Definition)) minBufferSize = sizeof(C##Name_##Definition);
350#define DECLARE_NODE_PAR(Name_, name_)
351#include "node_type.conf"
352#undef DECLARE_NODE
353#undef DECLARE_NODE_PAR
354
355     // Compute the buffer sizes needed to send the attributes and data corresponding to fields
356     std::map<int, StdSize> maxEventSize;
357     std::map<int, StdSize> bufferSize = getAttributesBufferSize(maxEventSize, contextClient, bufferForWriting);
358     std::map<int, StdSize> dataBufferSize = getDataBufferSize(maxEventSize, contextClient, bufferForWriting);
359
360     std::map<int, StdSize>::iterator it, ite = dataBufferSize.end();
361     for (it = dataBufferSize.begin(); it != ite; ++it)
362       if (it->second > bufferSize[it->first]) bufferSize[it->first] = it->second;
363
364     // Apply the buffer size factor, check that we are above the minimum buffer size and below the maximum size
365     ite = bufferSize.end();
366     for (it = bufferSize.begin(); it != ite; ++it)
367     {
368       it->second *= CXios::bufferSizeFactor;
369       if (it->second < minBufferSize) it->second = minBufferSize;
370       if (it->second > CXios::maxBufferSize) it->second = CXios::maxBufferSize;
371     }
372
373     // Leaders will have to send some control events so ensure there is some room for those in the buffers
374     if (contextClient->isServerLeader())
375     {
376       const std::list<int>& ranks = contextClient->getRanksServerLeader();
377       for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
378       {
379         if (!bufferSize.count(*itRank))
380         {
381           bufferSize[*itRank] = minBufferSize;
382           maxEventSize[*itRank] = minEventSize;
383         }
384       }
385     }
386     contextClient->setBufferSize(bufferSize, maxEventSize);
387   }
388   CATCH_DUMP_ATTR
389
390 /*!
391    * Compute the required buffer size to send the fields data.
392    * \param maxEventSize [in/out] the size of the bigger event for each connected server
393    * \param [in] contextClient
394    * \param [in] bufferForWriting True if buffers are used for sending data for writing
395      This flag is only true for client and server-1 for communication with server-2
396    */
397   std::map<int, StdSize> CContext::getDataBufferSize(std::map<int, StdSize>& maxEventSize,
398                                                      CContextClient* contextClient, bool bufferForWriting /*= "false"*/)
399   TRY
400   {
401     std::map<int, StdSize> dataSize;
402
403     // Find all reference domain and axis of all active fields
404     std::vector<CFile*>& fileList = bufferForWriting ? this->enabledWriteModeFiles : this->enabledReadModeFiles;
405     size_t numEnabledFiles = fileList.size();
406     for (size_t i = 0; i < numEnabledFiles; ++i)
407     {
408       CFile* file = fileList[i];
409       if (file->getContextClient() == contextClient)
410       {
411         std::vector<CField*> enabledFields = file->getEnabledFields();
412         size_t numEnabledFields = enabledFields.size();
413         for (size_t j = 0; j < numEnabledFields; ++j)
414         {
415           // const std::vector<std::map<int, StdSize> > mapSize = enabledFields[j]->getGridDataBufferSize(contextClient);
416           const std::map<int, StdSize> mapSize = enabledFields[j]->getGridDataBufferSize(contextClient,bufferForWriting);
417           std::map<int, StdSize>::const_iterator it = mapSize.begin(), itE = mapSize.end();
418           for (; it != itE; ++it)
419           {
420             // If dataSize[it->first] does not exist, it will be zero-initialized
421             // so we can use it safely without checking for its existance
422           if (CXios::isOptPerformance)
423               dataSize[it->first] += it->second;
424             else if (dataSize[it->first] < it->second)
425               dataSize[it->first] = it->second;
426
427           if (maxEventSize[it->first] < it->second)
428               maxEventSize[it->first] = it->second;
429           }
430         }
431       }
432     }
433     return dataSize;
434   }
435   CATCH_DUMP_ATTR
436
437/*!
438    * Compute the required buffer size to send the attributes (mostly those grid related).
439    * \param maxEventSize [in/out] the size of the bigger event for each connected server
440    * \param [in] contextClient
441    * \param [in] bufferForWriting True if buffers are used for sending data for writing
442      This flag is only true for client and server-1 for communication with server-2
443    */
444   std::map<int, StdSize> CContext::getAttributesBufferSize(std::map<int, StdSize>& maxEventSize,
445                                                           CContextClient* contextClient, bool bufferForWriting /*= "false"*/)
446   TRY
447   {
448   // As calendar attributes are sent even if there are no active files or fields, maps are initialized according the size of calendar attributes
449     std::map<int, StdSize> attributesSize = CCalendarWrapper::get(CCalendarWrapper::GetDefName())->getMinimumBufferSizeForAttributes(contextClient);
450     maxEventSize = CCalendarWrapper::get(CCalendarWrapper::GetDefName())->getMinimumBufferSizeForAttributes(contextClient);
451
452     std::vector<CFile*>& fileList = this->enabledFiles;
453     size_t numEnabledFiles = fileList.size();
454     for (size_t i = 0; i < numEnabledFiles; ++i)
455     {
456//         CFile* file = this->enabledWriteModeFiles[i];
457        CFile* file = fileList[i];
458        std::vector<CField*> enabledFields = file->getEnabledFields();
459        size_t numEnabledFields = enabledFields.size();
460        for (size_t j = 0; j < numEnabledFields; ++j)
461        {
462          const std::map<int, StdSize> mapSize = enabledFields[j]->getGridAttributesBufferSize(contextClient, bufferForWriting);
463          std::map<int, StdSize>::const_iterator it = mapSize.begin(), itE = mapSize.end();
464          for (; it != itE; ++it)
465          {
466         // If attributesSize[it->first] does not exist, it will be zero-initialized
467         // so we can use it safely without checking for its existence
468             if (attributesSize[it->first] < it->second)
469         attributesSize[it->first] = it->second;
470
471         if (maxEventSize[it->first] < it->second)
472         maxEventSize[it->first] = it->second;
473          }
474        }
475     }
476     return attributesSize;
477   }
478   CATCH_DUMP_ATTR
479
480
481
482   //! Verify whether a context is initialized
483   bool CContext::isInitialized(void)
484   TRY
485   {
486     return hasClient;
487   }
488   CATCH_DUMP_ATTR
489
490
491   void CContext::init(CServerContext* parentServerContext, MPI_Comm intraComm, int serviceType)
492   TRY
493   {
494     parentServerContext_ = parentServerContext ;
495     if (serviceType==CServicesManager::CLIENT) 
496       initClient(intraComm, serviceType) ;
497     else
498       initServer(intraComm, serviceType) ;
499    }
500    CATCH_DUMP_ATTR
501
502
503
504//! Initialize client side
505   void CContext::initClient(MPI_Comm intraComm, int serviceType)
506   TRY
507   {
508      intraComm_=intraComm ;
509      MPI_Comm_rank(intraComm_, &intraCommRank_) ;
510      MPI_Comm_size(intraComm_, &intraCommSize_) ;
511
512      serviceType_ = CServicesManager::CLIENT ;
513      if (serviceType_==CServicesManager::CLIENT)
514      {
515        hasClient=true ;
516        hasServer=false ;
517      }
518      contextId_ = getId() ;
519     
520      attached_mode=true ;
521      if (!CXios::isUsingServer()) attached_mode=false ;
522
523
524      string contextRegistryId=getId() ;
525      registryIn=new CRegistry(intraComm);
526      registryIn->setPath(contextRegistryId) ;
527     
528      int commRank ;
529      MPI_Comm_rank(intraComm_,&commRank) ;
530      if (commRank==0) registryIn->fromFile("xios_registry.bin") ;
531      registryIn->bcastRegistry() ;
532      registryOut=new CRegistry(intraComm_) ;
533      registryOut->setPath(contextRegistryId) ;
534     
535   }
536   CATCH_DUMP_ATTR
537
538   
539   void CContext::initServer(MPI_Comm intraComm, int serviceType)
540   TRY
541   {
542     hasServer=true;
543     intraComm_=intraComm ;
544     MPI_Comm_rank(intraComm_, &intraCommRank_) ;
545     MPI_Comm_size(intraComm_, &intraCommSize_) ;
546
547     serviceType_=serviceType ;
548
549     if (serviceType_==CServicesManager::GATHERER)
550     {
551       hasClient=true ;
552       hasServer=true ;
553     }
554     else if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)
555     {
556       hasClient=false ;
557       hasServer=true ;
558     }
559
560     CXios::getContextsManager()->getContextId(getId(), contextId_, intraComm) ;
561     
562     registryIn=new CRegistry(intraComm);
563     registryIn->setPath(contextId_) ;
564     
565     int commRank ;
566     MPI_Comm_rank(intraComm_,&commRank) ;
567     if (commRank==0) registryIn->fromFile("xios_registry.bin") ;
568   
569     registryIn->bcastRegistry() ;
570     registryOut=new CRegistry(intraComm) ;
571     registryOut->setPath(contextId_) ;
572
573   }
574   CATCH_DUMP_ATTR
575
576
577  void CContext::createClientInterComm(MPI_Comm interCommClient, MPI_Comm interCommServer) // for servers
578  TRY
579  {
580    MPI_Comm intraCommClient ;
581    MPI_Comm_dup(intraComm_, &intraCommClient);
582    comms.push_back(intraCommClient);
583    // attached_mode=parentServerContext_->isAttachedMode() ; //ym probably inherited from source context
584    server = new CContextServer(this,intraComm_, interCommServer); // check if we need to dupl. intraComm_ ?
585    client = new CContextClient(this,intraCommClient,interCommClient);
586    client->setAssociatedServer(server) ; 
587    server->setAssociatedClient(client) ; 
588
589  }
590  CATCH_DUMP_ATTR
591
592  void CContext::createServerInterComm(void) 
593  TRY
594  {
595   
596    MPI_Comm interCommClient, interCommServer ;
597
598    if (serviceType_ == CServicesManager::CLIENT)
599    {
600
601      int commRank ;
602      MPI_Comm_rank(intraComm_,&commRank) ;
603      if (commRank==0)
604      {
605        if (attached_mode) CXios::getContextsManager()->createServerContext(CClient::getPoolRessource()->getId(), CXios::defaultServerId, 0, getContextId()) ;
606        else if (CXios::usingServer2) CXios::getContextsManager()->createServerContext(CXios::defaultPoolId, CXios::defaultGathererId, 0, getContextId()) ;
607        else  CXios::getContextsManager()->createServerContext(CXios::defaultPoolId, CXios::defaultServerId, 0, getContextId()) ;
608      }
609
610      MPI_Comm interComm ;
611     
612      if (attached_mode)
613      {
614        parentServerContext_->createIntercomm(CClient::getPoolRessource()->getId(), CXios::defaultServerId, 0, getContextId(), intraComm_, 
615                                              interCommClient, interCommServer) ;
616        int type ; 
617        if (commRank==0) CXios::getServicesManager()->getServiceType(CClient::getPoolRessource()->getId(), CXios::defaultServerId, 0, type) ;
618        MPI_Bcast(&type,1,MPI_INT,0,intraComm_) ;
619        setCurrent(getId()) ; // getCurrent/setCurrent may be supress, it can cause a lot of trouble
620      }
621      else if (CXios::usingServer2)
622      { 
623//      CXios::getContextsManager()->createServerContextIntercomm(CXios::defaultPoolId, CXios::defaultGathererId, 0, getContextId(), intraComm_, interComm) ;
624        parentServerContext_->createIntercomm(CXios::defaultPoolId, CXios::defaultGathererId, 0, getContextId(), intraComm_,
625                                              interCommClient, interCommServer) ;
626        int type ; 
627        if (commRank==0) CXios::getServicesManager()->getServiceType(CXios::defaultPoolId, CXios::defaultGathererId, 0, type) ;
628        MPI_Bcast(&type,1,MPI_INT,0,intraComm_) ;
629      }
630      else
631      {
632        //CXios::getContextsManager()->createServerContextIntercomm(CXios::defaultPoolId, CXios::defaultServerId, 0, getContextId(), intraComm_, interComm) ;
633        parentServerContext_->createIntercomm(CXios::defaultPoolId, CXios::defaultServerId, 0, getContextId(), intraComm_,
634                                              interCommClient, interCommServer) ;
635        int type ; 
636        if (commRank==0) CXios::getServicesManager()->getServiceType(CXios::defaultPoolId, CXios::defaultServerId, 0, type) ;
637        MPI_Bcast(&type,1,MPI_INT,0,intraComm_) ;
638      }
639
640        // intraComm client is not duplicated. In all the code we use client->intraComm for MPI
641        // in future better to replace it by intracommuncator associated to the context
642   
643      MPI_Comm intraCommClient, intraCommServer ;
644      intraCommClient=intraComm_ ;
645      MPI_Comm_dup(intraComm_, &intraCommServer) ;
646      client = new CContextClient(this, intraCommClient, interCommClient);
647      server = new CContextServer(this, intraCommServer, interCommServer);
648      client->setAssociatedServer(server) ;
649      server->setAssociatedClient(client) ;
650    }
651   
652    if (serviceType_ == CServicesManager::GATHERER)
653    {
654      int commRank ;
655      MPI_Comm_rank(intraComm_,&commRank) ;
656     
657      int nbPartitions ;
658      if (commRank==0) 
659      { 
660        CXios::getServicesManager()->getServiceNbPartitions(CXios::defaultPoolId, CXios::defaultServerId, 0, nbPartitions) ;
661        for(int i=0 ; i<nbPartitions; i++)
662          CXios::getContextsManager()->createServerContext(CXios::defaultPoolId, CXios::defaultServerId, i, getContextId()) ;
663      }     
664      MPI_Bcast(&nbPartitions, 1, MPI_INT, 0, intraComm_) ;
665     
666      MPI_Comm interComm ;
667      for(int i=0 ; i<nbPartitions; i++)
668      {
669        parentServerContext_->createIntercomm(CXios::defaultPoolId, CXios::defaultServerId, i, getContextId(), intraComm_, interCommClient, interCommServer) ;
670        int type ; 
671        if (commRank==0) CXios::getServicesManager()->getServiceType(CXios::defaultPoolId, CXios::defaultServerId, 0, type) ;
672        MPI_Bcast(&type,1,MPI_INT,0,intraComm_) ;
673        primServerId_.push_back(CXios::getContextsManager()->getServerContextName(CXios::defaultPoolId, CXios::defaultServerId, i, type, getContextId())) ;
674
675        // intraComm client is not duplicated. In all the code we use client->intraComm for MPI
676        // in future better to replace it by intracommuncator associated to the context
677     
678        MPI_Comm intraCommClient, intraCommServer ;
679
680        intraCommClient=intraComm_ ;
681        MPI_Comm_dup(intraComm_, &intraCommServer) ;
682
683        CContextClient* client = new CContextClient(this, intraCommClient, interCommClient) ;
684        CContextServer* server = new CContextServer(this, intraCommServer, interCommServer) ;
685        client->setAssociatedServer(server) ;
686        server->setAssociatedClient(client) ;
687        clientPrimServer.push_back(client);
688        serverPrimServer.push_back(server); 
689
690     
691      }
692    }
693  }
694  CATCH_DUMP_ATTR
695
696   
697
698  bool CContext::eventLoop(bool enableEventsProcessing)
699  {
700    bool finished=true; 
701
702    if (client!=nullptr && !finalized) client->checkBuffers();
703   
704    for (int i = 0; i < clientPrimServer.size(); ++i)
705    {
706      if (!finalized) clientPrimServer[i]->checkBuffers();
707      if (!finalized) finished &= serverPrimServer[i]->eventLoop(enableEventsProcessing);
708    }
709
710    for (auto it=couplerClient_.begin(); it!=couplerClient_.end(); ++it)
711    {
712      if (!finalized) it->second->checkBuffers();
713    }
714
715    for (auto it=couplerServer_.begin(); it!=couplerServer_.end(); ++it)
716    {
717      if (!finalized) it->second->eventLoop(enableEventsProcessing);
718    }
719
720    if (server!=nullptr) if (!finalized) finished &= server->eventLoop(enableEventsProcessing);
721 
722    return finalized && finished ;
723  }
724
725  void CContext::addCouplingChanel(const std::string& context, bool out)
726  {
727     vector<string> vectStr=splitRegex(context,"::") ;
728     string poolId=vectStr[0] ;
729     string serviceId=poolId ;
730     string contextId=vectStr[1] ;
731
732     int contextLeader ;
733     int type = CServicesManager::CLIENT ;
734     string contextName=CXios::getContextsManager()->getServerContextName(poolId, serviceId, 0, type, contextId) ;
735     
736     if (couplerClient_.find(contextName)==couplerClient_.end())
737     {
738       bool ok=CXios::getContextsManager()->getContextLeader(contextName, contextLeader, getIntraComm()) ;
739     
740       MPI_Comm interComm, interCommClient, interCommServer  ;
741       MPI_Comm intraCommClient, intraCommServer ;
742
743       if (ok) MPI_Intercomm_create(getIntraComm(), 0, CXios::getXiosComm(), contextLeader, 0, &interComm) ;
744
745       MPI_Comm_dup(intraComm_, &intraCommClient) ;
746       MPI_Comm_dup(intraComm_, &intraCommServer) ;
747       if (out)
748       {
749         MPI_Comm_dup(interComm, &interCommClient) ;
750         MPI_Comm_dup(interComm, &interCommServer) ;
751         CContextClient* client = new CContextClient(this, intraCommClient, interCommClient);
752         CContextServer* server = new CContextServer(this, intraCommServer, interCommServer);
753         client->setAssociatedServer(server) ;
754         server->setAssociatedClient(client) ;
755       }
756       else
757       {
758          MPI_Comm_dup(interComm, &interCommServer) ;
759          MPI_Comm_dup(interComm, &interCommClient) ;
760          CContextServer* server = new CContextServer(this, intraCommServer, interCommServer);
761          CContextClient* client = new CContextClient(this, intraCommClient, interCommClient);
762          client->setAssociatedServer(server) ;
763          server->setAssociatedClient(client) ;
764       }
765       MPI_Comm_free(&interComm) ;
766
767
768      // for now, we don't now which beffer size must be used for client coupler
769      // It will be evaluated later. Fix a constant size for now...
770      // set to 10Mb for development
771       map<int,size_t> bufferSize, maxEventSize ;
772       for(int i=0;i<client->getRemoteSize();i++)
773       {
774         bufferSize[i]=10000000 ;
775         maxEventSize[i]=10000000 ;
776       }
777
778       client->setBufferSize(bufferSize, maxEventSize);   
779       
780       couplerClient_[contextName] = client ;
781       couplerServer_[contextName] = server ;
782     }
783  }
784 
785  void CContext::globalEventLoop(void)
786  {
787    CXios::getDaemonsManager()->eventLoop() ;
788    setCurrent(getId()) ;
789  }
790
791
792   void CContext::finalize(void)
793   TRY
794   {
795      registryOut->hierarchicalGatherRegistry() ;
796      if (server->intraCommRank==0) CXios::globalRegistry->mergeRegistry(*registryOut) ;
797
798      if (serviceType_==CServicesManager::CLIENT)
799      {
800        doPreTimestepOperationsForEnabledReadModeFiles(); // For now we only use server level 1 to read data
801
802        info(100)<<"DEBUG: context "<<getId()<<" Send client finalize"<<endl ;
803        client->finalize();
804        info(100)<<"DEBUG: context "<<getId()<<" Client finalize sent"<<endl ;
805        while (client->havePendingRequests()) client->checkBuffers();
806        info(100)<<"DEBUG: context "<<getId()<<" no pending request ok"<<endl ;
807        bool notifiedFinalized=false ;
808        do
809        {
810          notifiedFinalized=client->isNotifiedFinalized() ;
811        } while (!notifiedFinalized) ;
812        client->releaseBuffers();
813        info(100)<<"DEBUG: context "<<getId()<<" release client ok"<<endl ;
814      }
815      else if (serviceType_==CServicesManager::GATHERER)
816      {
817         for (int i = 0; i < clientPrimServer.size(); ++i)
818         {
819           clientPrimServer[i]->finalize();
820           bool bufferReleased;
821           do
822           {
823             clientPrimServer[i]->checkBuffers();
824             bufferReleased = !clientPrimServer[i]->havePendingRequests();
825           } while (!bufferReleased);
826           
827           bool notifiedFinalized=false ;
828           do
829           {
830             notifiedFinalized=clientPrimServer[i]->isNotifiedFinalized() ;
831           } while (!notifiedFinalized) ;
832           clientPrimServer[i]->releaseBuffers();
833         }
834         closeAllFile();
835
836      }
837      else if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)
838      {
839        closeAllFile();
840      }
841
842      freeComms() ;
843       
844      parentServerContext_->freeComm() ;
845      finalized = true;
846      info(20)<<"CContext: Context <"<<getId()<<"> is finalized."<<endl;
847   }
848   CATCH_DUMP_ATTR
849
850   //! Free internally allocated communicators
851   void CContext::freeComms(void)
852   TRY
853   {
854     for (std::list<MPI_Comm>::iterator it = comms.begin(); it != comms.end(); ++it)
855       MPI_Comm_free(&(*it));
856     comms.clear();
857   }
858   CATCH_DUMP_ATTR
859
860   //! Deallocate buffers allocated by clientContexts
861   void CContext::releaseClientBuffers(void)
862   TRY
863   {
864     client->releaseBuffers();
865     for (int i = 0; i < clientPrimServer.size(); ++i)
866       clientPrimServer[i]->releaseBuffers();
867   }
868   CATCH_DUMP_ATTR
869
870   void CContext::postProcessingGlobalAttributes()
871   TRY
872   {
873     if (allProcessed) return; 
874     
875    // create intercommunicator with servers.
876    // not sure it is the good place to be called here
877    createServerInterComm() ;
878
879
880     // After xml is parsed, there are some more works with post processing
881     postProcessing();
882
883     // Distribute files between secondary servers according to the data size
884     distributeFiles(this->enabledWriteModeFiles);
885
886     // Check grid and calculate its distribution
887     checkGridEnabledFields();
888
889     setClientServerBuffer(client, (serviceType_==CServicesManager::CLIENT) ) ;
890     for (int i = 0; i < clientPrimServer.size(); ++i)
891         setClientServerBuffer(clientPrimServer[i], true);
892
893   
894     if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
895     { 
896       if (serviceType_==CServicesManager::GATHERER)
897       { 
898         for (auto it=clientPrimServer.begin(); it!=clientPrimServer.end();++it) 
899         {
900           this->sendAllAttributesToServer(*it); // Send all attributes of current context to server
901           CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(*it); // Send all attributes of current calendar
902         }
903       }
904       else 
905       {
906         this->sendAllAttributesToServer(client);   // Send all attributes of current context to server
907         CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current calendar
908       }
909
910
911      // We have enough information to send to server
912      // First of all, send all enabled files
913      sendEnabledFiles(this->enabledWriteModeFiles);
914      // We only use server-level 1 (for now) to read data
915      if (serviceType_==CServicesManager::CLIENT)  sendEnabledFiles(this->enabledReadModeFiles);
916
917      // Then, send all enabled fields     
918      sendEnabledFieldsInFiles(this->enabledWriteModeFiles);
919     
920      if (serviceType_==CServicesManager::CLIENT) sendEnabledFieldsInFiles(this->enabledReadModeFiles);
921
922      // Then, check whether we have domain_ref, axis_ref or scalar_ref attached to the enabled fields
923      // If any, so send them to server
924       sendRefDomainsAxisScalars(this->enabledWriteModeFiles); 
925     
926      if (serviceType_==CServicesManager::CLIENT) sendRefDomainsAxisScalars(this->enabledReadModeFiles);       
927
928       // Check whether enabled fields have grid_ref, if any, send this info to server
929      sendRefGrid(this->enabledFiles);
930      // This code may be useful in the future when we want to seperate completely read and write
931      // sendRefGrid(this->enabledWriteModeFiles);
932      // if (!hasServer)
933      //   sendRefGrid(this->enabledReadModeFiles);
934     
935      // A grid of enabled fields composed of several components which must be checked then their
936      // checked attributes should be sent to server
937      sendGridComponentEnabledFieldsInFiles(this->enabledFiles); // This code can be seperated in two (one for reading, another for writing)
938
939       // We have a xml tree on the server side and now, it should be also processed
940      sendPostProcessing();
941       
942      // Finally, we send information of grid itself to server
943      sendGridEnabledFieldsInFiles(this->enabledWriteModeFiles);       
944     
945      if (serviceType_==CServicesManager::CLIENT) sendGridEnabledFieldsInFiles(this->enabledReadModeFiles);       
946     
947     }
948     allProcessed = true;
949   }
950   CATCH_DUMP_ATTR
951
952   void CContext::sendPostProcessingGlobalAttributes()
953   TRY
954   {
955
956    int nbSrvPools ;
957    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
958    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
959    else nbSrvPools = 0 ;
960    CContextClient* contextClientTmp ;
961
962    for (int i = 0; i < nbSrvPools; ++i)
963     {
964       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
965       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
966     
967       CEventClient event(getType(),EVENT_ID_POST_PROCESS_GLOBAL_ATTRIBUTES);
968
969       if (contextClientTmp->isServerLeader())
970       {
971         CMessage msg;
972         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
973         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
974           event.push(*itRank,1,msg);
975         contextClientTmp->sendEvent(event);
976       }
977       else contextClientTmp->sendEvent(event);
978     }
979   }
980   CATCH_DUMP_ATTR
981
982   void CContext::recvPostProcessingGlobalAttributes(CEventServer& event)
983   TRY
984   {
985      CBufferIn* buffer=event.subEvents.begin()->buffer;
986      getCurrent()->recvPostProcessingGlobalAttributes(*buffer);
987   }
988   CATCH
989
990   void CContext::recvPostProcessingGlobalAttributes(CBufferIn& buffer)
991   TRY
992   {     
993      postProcessingGlobalAttributes();
994   }
995   CATCH_DUMP_ATTR
996
997   /*!
998   \brief Close all the context defintion and do processing data
999      After everything is well defined on client side, they will be processed and sent to server
1000   From the version 2.0, sever and client work no more on the same database. Moreover, client(s) will send
1001   all necessary information to server, from which each server can build its own database.
1002   Because the role of server is to write out field data on a specific netcdf file,
1003   the only information that it needs is the enabled files
1004   and the active fields (fields will be written onto active files)
1005   */
1006  void CContext::closeDefinition(void)
1007   TRY
1008   {
1009     CTimer::get("Context : close definition").resume() ;
1010     
1011     // create intercommunicator with servers.
1012     // not sure it is the good place to be called here
1013     createServerInterComm() ;
1014
1015
1016     // After xml is parsed, there are some more works with post processing
1017//     postProcessing();
1018
1019    // Make sure the calendar was correctly created
1020    if (!calendar)
1021      ERROR("CContext::postProcessing()", << "A calendar must be defined for the context \"" << getId() << "!\"")
1022    else if (calendar->getTimeStep() == NoneDu)
1023      ERROR("CContext::postProcessing()", << "A timestep must be defined for the context \"" << getId() << "!\"")
1024    // Calendar first update to set the current date equals to the start date
1025    calendar->update(0);
1026
1027    // Résolution des héritages descendants (càd des héritages de groupes)
1028    // pour chacun des contextes.
1029    solveDescInheritance(true);
1030 
1031    // Solve inheritance for field to know if enabled or not.
1032    for (auto field : CField::getAll()) field->solveRefInheritance();
1033
1034    // Check if some axis, domains or grids are eligible to for compressed indexed output.
1035    // Warning: This must be done after solving the inheritance and before the rest of post-processing
1036    // --> later ????    checkAxisDomainsGridsEligibilityForCompressedOutput();     
1037
1038      // Check if some automatic time series should be generated
1039      // Warning: This must be done after solving the inheritance and before the rest of post-processing     
1040
1041    // The timeseries should only be prepared in client
1042    prepareTimeseries();
1043
1044    //Initialisation du vecteur 'enabledFiles' contenant la liste des fichiers à sortir.
1045    findEnabledFiles();
1046    findEnabledWriteModeFiles();
1047    findEnabledReadModeFiles();
1048    findEnabledCouplerIn();
1049    findEnabledCouplerOut();
1050    createCouplerInterCommunicator() ;
1051
1052    // Find all enabled fields of each file     
1053    vector<CField*>&& fileOutField = findAllEnabledFieldsInFileOut(this->enabledWriteModeFiles);
1054    vector<CField*>&& fileInField = findAllEnabledFieldsInFileIn(this->enabledReadModeFiles);
1055    vector<CField*>&& CouplerOutField = findAllEnabledFieldsCouplerOut(this->enabledCouplerOut);
1056    vector<CField*>&& CouplerInField = findAllEnabledFieldsCouplerIn(this->enabledCouplerIn);
1057    findFieldsWithReadAccess();
1058    vector<CField*>& fieldWithReadAccess = fieldsWithReadAccess_ ;
1059    vector<CField*> fieldModelIn ; // fields potentially from model
1060     
1061// find all field potentially at workflow end
1062    vector<CField*> endWorkflowFields ;
1063    endWorkflowFields.reserve(fileOutField.size()+CouplerOutField.size()+fieldWithReadAccess.size()) ;
1064    endWorkflowFields.insert(endWorkflowFields.end(),fileOutField.begin(), fileOutField.end()) ;
1065    endWorkflowFields.insert(endWorkflowFields.end(),CouplerOutField.begin(), CouplerOutField.end()) ;
1066    endWorkflowFields.insert(endWorkflowFields.end(),fieldWithReadAccess.begin(), fieldWithReadAccess.end()) ;
1067
1068    for(auto endWorkflowField : endWorkflowFields) endWorkflowField->buildWorkflowGraph(garbageCollector) ;
1069   
1070    // get all field coming potentially from model
1071    for (auto field : CField::getAll() ) if (field->getModelIn()) fieldModelIn.push_back(field) ;
1072
1073    // Distribute files between secondary servers according to the data size => assign a context to a file and then to fields
1074    if (serviceType_==CServicesManager::GATHERER) distributeFiles(this->enabledWriteModeFiles);
1075    else if (serviceType_==CServicesManager::CLIENT) for(auto file : this->enabledWriteModeFiles) file->setContextClient(client) ;
1076
1077   
1078    // workflow endpoint => sent to IO/SERVER
1079    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1080    {
1081      for(auto field : fileOutField) 
1082      {
1083        field->connectToFileServer(garbageCollector) ; // connect the field to server filter
1084        field->computeGridIndexToFileServer() ; // compute grid index for transfer to the server context
1085      }
1086      setClientServerBuffer(fileOutField, true) ; // set context
1087      for(auto field : fileOutField) field->sendFieldToFileServer() ;
1088    }
1089
1090
1091    // workflow startpoint => data from model
1092    if (serviceType_==CServicesManager::CLIENT)
1093    {
1094      for(auto field : fieldModelIn) 
1095      {
1096        field->connectToModelInput(garbageCollector) ; // connect the field to server filter
1097        // grid index will be computed on the fly
1098      }
1099    }
1100
1101 
1102
1103
1104    return ;
1105    // For now, only read files with client and only one level server
1106    // if (hasClient && !hasServer) findEnabledReadModeFiles();     
1107
1108    // Find all enabled fields of each file     
1109    findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
1110    findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
1111
1112    // For now, only read files with client and only one level server
1113    // if (hasClient && !hasServer)
1114    //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
1115
1116    if (serviceType_==CServicesManager::CLIENT)
1117    {
1118      initReadFiles();
1119      // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
1120      this->readAttributesOfEnabledFieldsInReadModeFiles();
1121    }
1122
1123    // Only search and rebuild all reference objects of enable fields, don't transform
1124    this->solveOnlyRefOfEnabledFields();
1125
1126    // Search and rebuild all reference object of enabled fields, and transform
1127    this->solveAllRefOfEnabledFieldsAndTransform();
1128
1129    // Find all fields with read access from the public API
1130    if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
1131    // and solve the all reference for them
1132    if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
1133
1134    isPostProcessed = true;
1135
1136
1137
1138    // Distribute files between secondary servers according to the data size
1139    distributeFiles(this->enabledWriteModeFiles);
1140
1141    // Check grid and calculate its distribution
1142    checkGridEnabledFields();
1143
1144    setClientServerBuffer(client, (serviceType_==CServicesManager::CLIENT) ) ;
1145    for (int i = 0; i < clientPrimServer.size(); ++i)
1146         setClientServerBuffer(clientPrimServer[i], true);
1147
1148   
1149    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1150    { 
1151      if (serviceType_==CServicesManager::GATHERER)
1152      { 
1153        for (auto it=clientPrimServer.begin(); it!=clientPrimServer.end();++it) 
1154        {
1155          this->sendAllAttributesToServer(*it); // Send all attributes of current context to server
1156          CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(*it); // Send all attributes of current calendar
1157        }
1158      }
1159      else 
1160      {
1161        this->sendAllAttributesToServer(client);   // Send all attributes of current context to server
1162        CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current calendar
1163      }
1164
1165      // We have enough information to send to server
1166      // First of all, send all enabled files
1167      sendEnabledFiles(this->enabledWriteModeFiles);
1168      // We only use server-level 1 (for now) to read data
1169      if (serviceType_==CServicesManager::CLIENT)  sendEnabledFiles(this->enabledReadModeFiles);
1170
1171      // Then, send all enabled fields     
1172      sendEnabledFieldsInFiles(this->enabledWriteModeFiles);
1173     
1174      if (serviceType_==CServicesManager::CLIENT) sendEnabledFieldsInFiles(this->enabledReadModeFiles);
1175
1176      // Then, check whether we have domain_ref, axis_ref or scalar_ref attached to the enabled fields
1177      // If any, so send them to server
1178      sendRefDomainsAxisScalars(this->enabledWriteModeFiles); 
1179     
1180      if (serviceType_==CServicesManager::CLIENT) sendRefDomainsAxisScalars(this->enabledReadModeFiles);       
1181
1182      // Check whether enabled fields have grid_ref, if any, send this info to server
1183      sendRefGrid(this->enabledFiles);
1184      // This code may be useful in the future when we want to seperate completely read and write
1185      // sendRefGrid(this->enabledWriteModeFiles);
1186      // if (!hasServer)
1187      //   sendRefGrid(this->enabledReadModeFiles);
1188     
1189      // A grid of enabled fields composed of several components which must be checked then their
1190      // checked attributes should be sent to server
1191      sendGridComponentEnabledFieldsInFiles(this->enabledFiles); // This code can be seperated in two (one for reading, another for writing)
1192
1193      // We have a xml tree on the server side and now, it should be also processed
1194      sendPostProcessing();
1195       
1196      // Finally, we send information of grid itself to server
1197      sendGridEnabledFieldsInFiles(this->enabledWriteModeFiles);       
1198     
1199      if (serviceType_==CServicesManager::CLIENT) sendGridEnabledFieldsInFiles(this->enabledReadModeFiles);       
1200    }
1201    allProcessed = true;
1202
1203
1204    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1205
1206    // There are some processings that should be done after all of above. For example: check mask or index
1207    this->buildFilterGraphOfEnabledFields();
1208   
1209     if (serviceType_==CServicesManager::CLIENT)
1210    {
1211      buildFilterGraphOfFieldsWithReadAccess();
1212      postProcessFilterGraph(); // For coupling in, modify this later
1213    }
1214   
1215    checkGridEnabledFields();   
1216
1217    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1218    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1219
1220    // Nettoyage de l'arborescence
1221    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1222
1223    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1224    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1225   
1226    CTimer::get("Context : close definition").suspend() ;
1227  }
1228  CATCH_DUMP_ATTR
1229
1230 /*!
1231  * Send context attribute and calendar to file server, it must be done once by context file server
1232  * \param[in] client : context client to send   
1233  */ 
1234  void CContext::sendContextToFileServer(CContextClient* client)
1235  {
1236    if (sendToFileServer_done_.count(client)!=0) return ;
1237    else sendToFileServer_done_.insert(client) ;
1238   
1239    this->sendAllAttributesToServer(client); // Send all attributes of current context to server
1240    CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current cale
1241  }
1242
1243  // ym obsolete now to be removed
1244   void CContext::closeDefinition_old(void)
1245   TRY
1246   {
1247     CTimer::get("Context : close definition").resume() ;
1248   
1249    //
1250    postProcessingGlobalAttributes();
1251
1252    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1253
1254    // There are some processings that should be done after all of above. For example: check mask or index
1255    this->buildFilterGraphOfEnabledFields();
1256   
1257     if (serviceType_==CServicesManager::CLIENT)
1258    {
1259      buildFilterGraphOfFieldsWithReadAccess();
1260      postProcessFilterGraph(); // For coupling in, modify this later
1261    }
1262   
1263    checkGridEnabledFields();   
1264
1265    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1266    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1267
1268    // Nettoyage de l'arborescence
1269    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1270
1271    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1272    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1273   
1274    CTimer::get("Context : close definition").suspend() ;
1275   }
1276   CATCH_DUMP_ATTR
1277
1278   vector<CField*> CContext::findAllEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1279   TRY
1280   {
1281     vector<CField*> fields ;
1282     for (unsigned int i = 0; i < activeFiles.size(); i++)
1283     {
1284        const vector<CField*>&& field=activeFiles[i]->getEnabledFields() ;
1285        fields.insert(fields.end(),field.begin(),field.end());
1286     }
1287     return fields ;
1288   }
1289   CATCH_DUMP_ATTR
1290
1291   vector<CField*> CContext::findAllEnabledFieldsInFileOut(const std::vector<CFile*>& activeFiles)
1292   TRY
1293   {
1294     vector<CField*> fields ;
1295     for(auto file : activeFiles)
1296     {
1297        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1298        for(auto field : fieldList) field->setFileOut(file) ;
1299        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1300     }
1301     return fields ;
1302   }
1303   CATCH_DUMP_ATTR
1304
1305   vector<CField*> CContext::findAllEnabledFieldsInFileIn(const std::vector<CFile*>& activeFiles)
1306   TRY
1307   {
1308     vector<CField*> fields ;
1309     for(auto file : activeFiles)
1310     {
1311        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1312        for(auto field : fieldList) field->setFileIn(file) ;
1313        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1314     }
1315     return fields ;
1316   }
1317   CATCH_DUMP_ATTR
1318
1319   vector<CField*> CContext::findAllEnabledFieldsCouplerOut(const std::vector<CCouplerOut*>& activeCouplerOut)
1320   TRY
1321   {
1322     vector<CField*> fields ;
1323     for (auto couplerOut :activeCouplerOut)
1324     {
1325        const vector<CField*>&& fieldList=couplerOut->getEnabledFields() ;
1326        for(auto field : fieldList) field->setCouplerOut(couplerOut) ;
1327        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1328     }
1329     return fields ;
1330   }
1331   CATCH_DUMP_ATTR
1332
1333   vector<CField*> CContext::findAllEnabledFieldsCouplerIn(const std::vector<CCouplerIn*>& activeCouplerIn)
1334   TRY
1335   {
1336     vector<CField*> fields ;
1337     for (auto couplerIn :activeCouplerIn)
1338     {
1339        const vector<CField*>&& fieldList=couplerIn->getEnabledFields() ;
1340        for(auto field : fieldList) field->setCouplerIn(couplerIn) ;
1341        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1342     }
1343     return fields ;
1344   }
1345   CATCH_DUMP_ATTR
1346
1347
1348
1349   void CContext::readAttributesOfEnabledFieldsInReadModeFiles()
1350   TRY
1351   {
1352      for (unsigned int i = 0; i < this->enabledReadModeFiles.size(); ++i)
1353        (void)this->enabledReadModeFiles[i]->readAttributesOfEnabledFieldsInReadMode();
1354   }
1355   CATCH_DUMP_ATTR
1356
1357   void CContext::sendGridComponentEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1358   TRY
1359   {
1360     int size = activeFiles.size();
1361     for (int i = 0; i < size; ++i)
1362     {       
1363       activeFiles[i]->sendGridComponentOfEnabledFields();
1364     }
1365   }
1366   CATCH_DUMP_ATTR
1367
1368   /*!
1369      Send active (enabled) fields in file from a client to others
1370      \param [in] activeFiles files contains enabled fields to send
1371   */
1372   void CContext::sendGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1373   TRY
1374   {
1375     int size = activeFiles.size();
1376     for (int i = 0; i < size; ++i)
1377     {       
1378       activeFiles[i]->sendGridOfEnabledFields();
1379     }
1380   }
1381   CATCH_DUMP_ATTR
1382
1383   void CContext::checkGridEnabledFields()
1384   TRY
1385   {
1386     int size = enabledFiles.size();
1387     for (int i = 0; i < size; ++i)
1388     {
1389       enabledFiles[i]->checkGridOfEnabledFields();       
1390     }
1391
1392     size = enabledCouplerOut.size();
1393     for (int i = 0; i < size; ++i)
1394     {
1395       enabledCouplerOut[i]->checkGridOfEnabledFields();       
1396     }
1397   }
1398   CATCH_DUMP_ATTR
1399
1400   /*!
1401      Check grid of active (enabled) fields in file
1402      \param [in] activeFiles files contains enabled fields whose grid needs checking
1403   */
1404   void CContext::checkGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1405   TRY
1406   {
1407     int size = activeFiles.size();
1408     for (int i = 0; i < size; ++i)
1409     {
1410       activeFiles[i]->checkGridOfEnabledFields();       
1411     }
1412   }
1413   CATCH_DUMP_ATTR
1414
1415    /*!
1416      Go up the hierachical tree via field_ref and do check of attributes of fields
1417      This can be done in a client then all computed information will be sent from this client to others
1418      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1419   */
1420   void CContext::solveOnlyRefOfEnabledFields(void)
1421   TRY
1422   {
1423     int size = this->enabledFiles.size();
1424     for (int i = 0; i < size; ++i)
1425     {
1426       this->enabledFiles[i]->solveOnlyRefOfEnabledFields();
1427     }
1428
1429     for (int i = 0; i < size; ++i)
1430     {
1431       this->enabledFiles[i]->generateNewTransformationGridDest();
1432     }
1433
1434     size = this->enabledCouplerOut.size();
1435     for (int i = 0; i < size; ++i)
1436     {
1437       this->enabledCouplerOut[i]->solveOnlyRefOfEnabledFields();
1438     }
1439
1440     for (int i = 0; i < size; ++i)
1441     {
1442       this->enabledCouplerOut[i]->generateNewTransformationGridDest();
1443     }
1444   }
1445   CATCH_DUMP_ATTR
1446
1447    /*!
1448      Go up the hierachical tree via field_ref and do check of attributes of fields.
1449      The transformation can be done in this step.
1450      All computed information will be sent from this client to others.
1451      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1452   */
1453   void CContext::solveAllRefOfEnabledFieldsAndTransform(void)
1454   TRY
1455   {
1456     int size = this->enabledFiles.size();
1457     for (int i = 0; i < size; ++i)
1458     {
1459       this->enabledFiles[i]->solveAllRefOfEnabledFieldsAndTransform();
1460     }
1461
1462     size = this->enabledCouplerOut.size();
1463     for (int i = 0; i < size; ++i)
1464     {
1465       this->enabledCouplerOut[i]->solveAllRefOfEnabledFieldsAndTransform();
1466     }
1467
1468   }
1469   CATCH_DUMP_ATTR
1470
1471   void CContext::buildFilterGraphOfEnabledFields()
1472   TRY
1473   {
1474     int size = this->enabledFiles.size();
1475     for (int i = 0; i < size; ++i)
1476     {
1477       this->enabledFiles[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1478     }
1479
1480     size = this->enabledCouplerOut.size();
1481     for (int i = 0; i < size; ++i)
1482     {
1483       this->enabledCouplerOut[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1484     }
1485   }
1486   CATCH_DUMP_ATTR
1487
1488   void CContext::postProcessFilterGraph()
1489   TRY
1490   {
1491     int size = enabledFiles.size();
1492     for (int i = 0; i < size; ++i)
1493     {
1494        enabledFiles[i]->postProcessFilterGraph();
1495     }
1496   }
1497   CATCH_DUMP_ATTR
1498
1499   void CContext::startPrefetchingOfEnabledReadModeFiles()
1500   TRY
1501   {
1502     int size = enabledReadModeFiles.size();
1503     for (int i = 0; i < size; ++i)
1504     {
1505        enabledReadModeFiles[i]->prefetchEnabledReadModeFields();
1506     }
1507   }
1508   CATCH_DUMP_ATTR
1509
1510   void CContext::doPreTimestepOperationsForEnabledReadModeFiles()
1511   TRY
1512   {
1513     int size = enabledReadModeFiles.size();
1514     for (int i = 0; i < size; ++i)
1515     {
1516        enabledReadModeFiles[i]->doPreTimestepOperationsForEnabledReadModeFields();
1517     }
1518   }
1519   CATCH_DUMP_ATTR
1520
1521   void CContext::doPostTimestepOperationsForEnabledReadModeFiles()
1522   TRY
1523   {
1524     int size = enabledReadModeFiles.size();
1525     for (int i = 0; i < size; ++i)
1526     {
1527        enabledReadModeFiles[i]->doPostTimestepOperationsForEnabledReadModeFields();
1528     }
1529   }
1530   CATCH_DUMP_ATTR
1531
1532  void CContext::findFieldsWithReadAccess(void)
1533  TRY
1534  {
1535    fieldsWithReadAccess_.clear();
1536    const vector<CField*> allFields = CField::getAll();
1537    for (size_t i = 0; i < allFields.size(); ++i)
1538    {
1539      CField* field = allFields[i];
1540      if (!field->read_access.isEmpty() && field->read_access && (field->enabled.isEmpty() || field->enabled))
1541      {
1542        fieldsWithReadAccess_.push_back(field);
1543        field->setModelOut() ;
1544      }
1545    }
1546  }
1547  CATCH_DUMP_ATTR
1548
1549  void CContext::solveAllRefOfFieldsWithReadAccess()
1550  TRY
1551  {
1552    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1553      fieldsWithReadAccess_[i]->solveAllReferenceEnabledField(false);
1554  }
1555  CATCH_DUMP_ATTR
1556
1557  void CContext::buildFilterGraphOfFieldsWithReadAccess()
1558  TRY
1559  {
1560    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1561      fieldsWithReadAccess_[i]->buildFilterGraph(garbageCollector, true);
1562  }
1563  CATCH_DUMP_ATTR
1564
1565   void CContext::solveAllInheritance(bool apply)
1566   TRY
1567   {
1568     // Résolution des héritages descendants (càd des héritages de groupes)
1569     // pour chacun des contextes.
1570      solveDescInheritance(apply);
1571
1572     // Résolution des héritages par référence au niveau des fichiers.
1573      const vector<CFile*> allFiles=CFile::getAll();
1574      const vector<CCouplerIn*> allCouplerIn=CCouplerIn::getAll();
1575      const vector<CCouplerOut*> allCouplerOut=CCouplerOut::getAll();
1576      const vector<CGrid*> allGrids= CGrid::getAll();
1577
1578      if (serviceType_==CServicesManager::CLIENT)
1579      {
1580        for (unsigned int i = 0; i < allFiles.size(); i++)
1581          allFiles[i]->solveFieldRefInheritance(apply);
1582
1583        for (unsigned int i = 0; i < allCouplerIn.size(); i++)
1584          allCouplerIn[i]->solveFieldRefInheritance(apply);
1585
1586        for (unsigned int i = 0; i < allCouplerOut.size(); i++)
1587          allCouplerOut[i]->solveFieldRefInheritance(apply);
1588      }
1589
1590      unsigned int vecSize = allGrids.size();
1591      unsigned int i = 0;
1592      for (i = 0; i < vecSize; ++i)
1593        allGrids[i]->solveElementsRefInheritance(apply);
1594
1595   }
1596  CATCH_DUMP_ATTR
1597
1598   void CContext::findEnabledFiles(void)
1599   TRY
1600   {
1601      const std::vector<CFile*> allFiles = CFile::getAll();
1602      const CDate& initDate = calendar->getInitDate();
1603
1604      for (unsigned int i = 0; i < allFiles.size(); i++)
1605         if (!allFiles[i]->enabled.isEmpty()) // Si l'attribut 'enabled' est défini.
1606         {
1607            if (allFiles[i]->enabled.getValue()) // Si l'attribut 'enabled' est fixé à vrai.
1608            {
1609              if (allFiles[i]->output_freq.isEmpty())
1610              {
1611                 ERROR("CContext::findEnabledFiles()",
1612                     << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1613                     <<" \".")
1614              }
1615              if ((initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1616              {
1617                error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1618                    << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1619                    <<"\" is less than the time step. File will not be written."<<endl;
1620              }
1621              else
1622               enabledFiles.push_back(allFiles[i]);
1623            }
1624         }
1625         else
1626         {
1627           if (allFiles[i]->output_freq.isEmpty())
1628           {
1629              ERROR("CContext::findEnabledFiles()",
1630                  << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1631                  <<" \".")
1632           }
1633           if ( (initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1634           {
1635             error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1636                 << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1637                 <<"\" is less than the time step. File will not be written."<<endl;
1638           }
1639           else
1640             enabledFiles.push_back(allFiles[i]); // otherwise true by default
1641         }
1642
1643      if (enabledFiles.size() == 0)
1644         DEBUG(<<"Aucun fichier ne va être sorti dans le contexte nommé \""
1645               << getId() << "\" !");
1646
1647   }
1648   CATCH_DUMP_ATTR
1649
1650   void CContext::findEnabledCouplerIn(void)
1651   TRY
1652   {
1653      const std::vector<CCouplerIn*> allCouplerIn = CCouplerIn::getAll();
1654      bool enabled ;
1655      for (size_t i = 0; i < allCouplerIn.size(); i++)
1656      {
1657        if (allCouplerIn[i]->enabled.isEmpty()) enabled=true ;
1658        else enabled=allCouplerIn[i]->enabled ;
1659        if (enabled) enabledCouplerIn.push_back(allCouplerIn[i]) ;
1660      }
1661   }
1662   CATCH_DUMP_ATTR
1663
1664   void CContext::findEnabledCouplerOut(void)
1665   TRY
1666   {
1667      const std::vector<CCouplerOut*> allCouplerOut = CCouplerOut::getAll();
1668      bool enabled ;
1669      for (size_t i = 0; i < allCouplerOut.size(); i++)
1670      {
1671        if (allCouplerOut[i]->enabled.isEmpty()) enabled=true ;
1672        else enabled=allCouplerOut[i]->enabled ;
1673        if (enabled) enabledCouplerOut.push_back(allCouplerOut[i]) ;
1674      }
1675   }
1676   CATCH_DUMP_ATTR
1677
1678
1679
1680
1681   void CContext::distributeFiles(const vector<CFile*>& files)
1682   TRY
1683   {
1684     bool distFileMemory=false ;
1685     distFileMemory=CXios::getin<bool>("server2_dist_file_memory", distFileMemory);
1686
1687     if (distFileMemory) distributeFileOverMemoryBandwith(files) ;
1688     else distributeFileOverBandwith(files) ;
1689   }
1690   CATCH_DUMP_ATTR
1691
1692   void CContext::distributeFileOverBandwith(const vector<CFile*>& files)
1693   TRY
1694   {
1695     double eps=std::numeric_limits<double>::epsilon()*10 ;
1696     
1697     std::ofstream ofs(("distribute_file_"+getId()+".dat").c_str(), std::ofstream::out);
1698     int nbPools = clientPrimServer.size();
1699
1700     // (1) Find all enabled files in write mode
1701     // for (int i = 0; i < this->enabledFiles.size(); ++i)
1702     // {
1703     //   if (enabledFiles[i]->mode.isEmpty() || (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
1704     //    enabledWriteModeFiles.push_back(enabledFiles[i]);
1705     // }
1706
1707     // (2) Estimate the data volume for each file
1708     int size = files.size();
1709     std::vector<std::pair<double, CFile*> > dataSizeMap;
1710     double dataPerPool = 0;
1711     int nfield=0 ;
1712     ofs<<size<<endl ;
1713     for (size_t i = 0; i < size; ++i)
1714     {
1715       CFile* file = files[i];
1716       ofs<<file->getId()<<endl ;
1717       StdSize dataSize=0;
1718       std::vector<CField*> enabledFields = file->getEnabledFields();
1719       size_t numEnabledFields = enabledFields.size();
1720       ofs<<numEnabledFields<<endl ;
1721       for (size_t j = 0; j < numEnabledFields; ++j)
1722       {
1723         dataSize += enabledFields[j]->getGlobalWrittenSize() ;
1724         ofs<<enabledFields[j]->getGrid()->getId()<<endl ;
1725         ofs<<enabledFields[j]->getGlobalWrittenSize()<<endl ;
1726       }
1727       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
1728       double dataSizeSec= dataSize/ outFreqSec;
1729       ofs<<dataSizeSec<<endl ;
1730       nfield++ ;
1731// add epsilon*nField to dataSizeSec in order to  preserve reproductive ordering when sorting
1732       dataSizeMap.push_back(make_pair(dataSizeSec + dataSizeSec * eps * nfield , file));
1733       dataPerPool += dataSizeSec;
1734     }
1735     dataPerPool /= nbPools;
1736     std::sort(dataSizeMap.begin(), dataSizeMap.end());
1737
1738     // (3) Assign contextClient to each enabled file
1739
1740     std::multimap<double,int> poolDataSize ;
1741// multimap is not garanty to preserve stable sorting in c++98 but it seems it does for c++11
1742
1743     int j;
1744     double dataSize ;
1745     for (j = 0 ; j < nbPools ; ++j) poolDataSize.insert(std::pair<double,int>(0.,j)) ; 
1746             
1747     for (int i = dataSizeMap.size()-1; i >= 0; --i)
1748     {
1749       dataSize=(*poolDataSize.begin()).first ;
1750       j=(*poolDataSize.begin()).second ;
1751       dataSizeMap[i].second->setContextClient(clientPrimServer[j]);
1752       dataSize+=dataSizeMap[i].first;
1753       poolDataSize.erase(poolDataSize.begin()) ;
1754       poolDataSize.insert(std::pair<double,int>(dataSize,j)) ; 
1755     }
1756
1757     for (std::multimap<double,int>:: iterator it=poolDataSize.begin() ; it!=poolDataSize.end(); ++it) info(30)<<"Load Balancing for servers (perfect=1) : "<<it->second<<" :  ratio "<<it->first*1./dataPerPool<<endl ;
1758   }
1759   CATCH_DUMP_ATTR
1760
1761   void CContext::distributeFileOverMemoryBandwith(const vector<CFile*>& filesList)
1762   TRY
1763   {
1764     int nbPools = clientPrimServer.size();
1765     double ratio=0.5 ;
1766     ratio=CXios::getin<double>("server2_dist_file_memory_ratio", ratio);
1767
1768     int nFiles = filesList.size();
1769     vector<SDistFile> files(nFiles);
1770     vector<SDistGrid> grids;
1771     map<string,int> gridMap ;
1772     string gridId; 
1773     int gridIndex=0 ;
1774
1775     for (size_t i = 0; i < nFiles; ++i)
1776     {
1777       StdSize dataSize=0;
1778       CFile* file = filesList[i];
1779       std::vector<CField*> enabledFields = file->getEnabledFields();
1780       size_t numEnabledFields = enabledFields.size();
1781
1782       files[i].id_=file->getId() ;
1783       files[i].nbGrids_=numEnabledFields;
1784       files[i].assignedGrid_ = new int[files[i].nbGrids_] ;
1785         
1786       for (size_t j = 0; j < numEnabledFields; ++j)
1787       {
1788         gridId=enabledFields[j]->getGrid()->getId() ;
1789         if (gridMap.find(gridId)==gridMap.end())
1790         {
1791            gridMap[gridId]=gridIndex  ;
1792            SDistGrid newGrid; 
1793            grids.push_back(newGrid) ;
1794            gridIndex++ ;
1795         }
1796         files[i].assignedGrid_[j]=gridMap[gridId] ;
1797         grids[files[i].assignedGrid_[j]].size_=enabledFields[j]->getGlobalWrittenSize() ;
1798         dataSize += enabledFields[j]->getGlobalWrittenSize() ; // usefull
1799       }
1800       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
1801       files[i].bandwith_= dataSize/ outFreqSec ;
1802     }
1803
1804     double bandwith=0 ;
1805     double memory=0 ;
1806   
1807     for(int i=0; i<nFiles; i++)  bandwith+=files[i].bandwith_ ;
1808     for(int i=0; i<nFiles; i++)  files[i].bandwith_ = files[i].bandwith_/bandwith * ratio ;
1809
1810     for(int i=0; i<grids.size(); i++)  memory+=grids[i].size_ ;
1811     for(int i=0; i<grids.size(); i++)  grids[i].size_ = grids[i].size_ / memory * (1.0-ratio) ;
1812       
1813     distributeFileOverServer2(nbPools, grids.size(), &grids[0], nFiles, &files[0]) ;
1814
1815     vector<double> memorySize(nbPools,0.) ;
1816     vector< set<int> > serverGrids(nbPools) ;
1817     vector<double> bandwithSize(nbPools,0.) ;
1818       
1819     for (size_t i = 0; i < nFiles; ++i)
1820     {
1821       bandwithSize[files[i].assignedServer_] += files[i].bandwith_* bandwith /ratio ;
1822       for(int j=0 ; j<files[i].nbGrids_;j++)
1823       {
1824         if (serverGrids[files[i].assignedServer_].find(files[i].assignedGrid_[j]) == serverGrids[files[i].assignedServer_].end())
1825         {
1826           memorySize[files[i].assignedServer_]+= grids[files[i].assignedGrid_[j]].size_ * memory / (1.0-ratio);
1827           serverGrids[files[i].assignedServer_].insert(files[i].assignedGrid_[j]) ;
1828         }
1829       }
1830       filesList[i]->setContextClient(clientPrimServer[files[i].assignedServer_]) ;
1831       delete [] files[i].assignedGrid_ ;
1832     }
1833
1834     for (int i = 0; i < nbPools; ++i) info(100)<<"Pool server level2 "<<i<<"   assigned file bandwith "<<bandwithSize[i]*86400.*4./1024/1024.<<" Mb / days"<<endl ;
1835     for (int i = 0; i < nbPools; ++i) info(100)<<"Pool server level2 "<<i<<"   assigned grid memory "<<memorySize[i]*100/1024./1024.<<" Mb"<<endl ;
1836
1837   }
1838   CATCH_DUMP_ATTR
1839
1840   /*!
1841      Find all files in write mode
1842   */
1843   void CContext::findEnabledWriteModeFiles(void)
1844   TRY
1845   {
1846     int size = this->enabledFiles.size();
1847     for (int i = 0; i < size; ++i)
1848     {
1849       if (enabledFiles[i]->mode.isEmpty() || 
1850          (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
1851        enabledWriteModeFiles.push_back(enabledFiles[i]);
1852     }
1853   }
1854   CATCH_DUMP_ATTR
1855
1856   /*!
1857      Find all files in read mode
1858   */
1859   void CContext::findEnabledReadModeFiles(void)
1860   TRY
1861   {
1862     int size = this->enabledFiles.size();
1863     for (int i = 0; i < size; ++i)
1864     {
1865       if (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::read)
1866        enabledReadModeFiles.push_back(enabledFiles[i]);
1867     }
1868   }
1869   CATCH_DUMP_ATTR
1870
1871   void CContext::closeAllFile(void)
1872   TRY
1873   {
1874     std::vector<CFile*>::const_iterator
1875            it = this->enabledFiles.begin(), end = this->enabledFiles.end();
1876
1877     for (; it != end; it++)
1878     {
1879       info(30)<<"Closing File : "<<(*it)->getId()<<endl;
1880       (*it)->close();
1881     }
1882   }
1883   CATCH_DUMP_ATTR
1884
1885   /*!
1886   \brief Dispatch event received from client
1887      Whenever a message is received in buffer of server, it will be processed depending on
1888   its event type. A new event type should be added in the switch list to make sure
1889   it processed on server side.
1890   \param [in] event: Received message
1891   */
1892   bool CContext::dispatchEvent(CEventServer& event)
1893   TRY
1894   {
1895
1896      if (SuperClass::dispatchEvent(event)) return true;
1897      else
1898      {
1899        switch(event.type)
1900        {
1901           case EVENT_ID_CLOSE_DEFINITION :
1902             recvCloseDefinition(event);
1903             return true;
1904             break;
1905           case EVENT_ID_UPDATE_CALENDAR:
1906             recvUpdateCalendar(event);
1907             return true;
1908             break;
1909           case EVENT_ID_CREATE_FILE_HEADER :
1910             recvCreateFileHeader(event);
1911             return true;
1912             break;
1913           case EVENT_ID_POST_PROCESS:
1914             recvPostProcessing(event);
1915             return true;
1916            case EVENT_ID_SEND_REGISTRY:
1917             recvRegistry(event);
1918             return true;
1919             break;
1920            case EVENT_ID_POST_PROCESS_GLOBAL_ATTRIBUTES:
1921             recvPostProcessingGlobalAttributes(event);
1922             return true;
1923             break;
1924            case EVENT_ID_PROCESS_GRID_ENABLED_FIELDS:
1925             recvProcessingGridOfEnabledFields(event);
1926             return true;
1927             break;
1928           default :
1929             ERROR("bool CContext::dispatchEvent(CEventServer& event)",
1930                    <<"Unknown Event");
1931           return false;
1932         }
1933      }
1934   }
1935   CATCH
1936
1937   //! Client side: Send a message to server to make it close
1938   void CContext::sendCloseDefinition(void)
1939   TRY
1940   {
1941    int nbSrvPools ;
1942    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
1943    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
1944    else nbSrvPools = 0 ;
1945    CContextClient* contextClientTmp ;
1946
1947    for (int i = 0; i < nbSrvPools; ++i)
1948     {
1949       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
1950       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
1951       CEventClient event(getType(),EVENT_ID_CLOSE_DEFINITION);
1952       if (contextClientTmp->isServerLeader())
1953       {
1954         CMessage msg;
1955         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
1956         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1957           event.push(*itRank,1,msg);
1958         contextClientTmp->sendEvent(event);
1959       }
1960       else contextClientTmp->sendEvent(event);
1961     }
1962   }
1963   CATCH_DUMP_ATTR
1964
1965   //! Server side: Receive a message of client announcing a context close
1966   void CContext::recvCloseDefinition(CEventServer& event)
1967   TRY
1968   {
1969      CBufferIn* buffer=event.subEvents.begin()->buffer;
1970      getCurrent()->closeDefinition();
1971   }
1972   CATCH
1973
1974   //! Client side: Send a message to update calendar in each time step
1975   void CContext::sendUpdateCalendar(int step)
1976   TRY
1977   {
1978    int nbSrvPools ;
1979    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
1980    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
1981    else nbSrvPools = 0 ;
1982    CContextClient* contextClientTmp ;
1983
1984    for (int i = 0; i < nbSrvPools; ++i)
1985    {
1986       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
1987       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
1988       CEventClient event(getType(),EVENT_ID_UPDATE_CALENDAR);
1989
1990         if (contextClientTmp->isServerLeader())
1991         {
1992           CMessage msg;
1993           msg<<step;
1994           const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
1995           for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1996             event.push(*itRank,1,msg);
1997           contextClientTmp->sendEvent(event);
1998         }
1999         else contextClientTmp->sendEvent(event);
2000     }
2001   }
2002   CATCH_DUMP_ATTR
2003
2004   //! Server side: Receive a message of client annoucing calendar update
2005   void CContext::recvUpdateCalendar(CEventServer& event)
2006   TRY
2007   {
2008      CBufferIn* buffer=event.subEvents.begin()->buffer;
2009      getCurrent()->recvUpdateCalendar(*buffer);
2010   }
2011   CATCH
2012
2013   //! Server side: Receive a message of client annoucing calendar update
2014   void CContext::recvUpdateCalendar(CBufferIn& buffer)
2015   TRY
2016   {
2017      int step;
2018      buffer>>step;
2019      updateCalendar(step);
2020      if (serviceType_==CServicesManager::GATHERER)
2021      {       
2022        sendUpdateCalendar(step);
2023      }
2024   }
2025   CATCH_DUMP_ATTR
2026
2027   //! Client side: Send a message to create header part of netcdf file
2028   void CContext::sendCreateFileHeader(void)
2029   TRY
2030   {
2031     int nbSrvPools ;
2032     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2033     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2034     else nbSrvPools = 0 ;
2035     CContextClient* contextClientTmp ;
2036
2037     for (int i = 0; i < nbSrvPools; ++i)
2038     {
2039       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2040       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2041       CEventClient event(getType(),EVENT_ID_CREATE_FILE_HEADER);
2042
2043       if (contextClientTmp->isServerLeader())
2044       {
2045         CMessage msg;
2046         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2047         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2048           event.push(*itRank,1,msg) ;
2049         contextClientTmp->sendEvent(event);
2050       }
2051       else contextClientTmp->sendEvent(event);
2052     }
2053   }
2054   CATCH_DUMP_ATTR
2055
2056   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2057   void CContext::recvCreateFileHeader(CEventServer& event)
2058   TRY
2059   {
2060      CBufferIn* buffer=event.subEvents.begin()->buffer;
2061      getCurrent()->recvCreateFileHeader(*buffer);
2062   }
2063   CATCH
2064
2065   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2066   void CContext::recvCreateFileHeader(CBufferIn& buffer)
2067   TRY
2068   {
2069      if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER) 
2070        createFileHeader();
2071   }
2072   CATCH_DUMP_ATTR
2073
2074   //! Client side: Send a message to do some post processing on server
2075   void CContext::sendProcessingGridOfEnabledFields()
2076   TRY
2077   {
2078     int nbSrvPools ;
2079     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2080     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2081     else nbSrvPools = 0 ;
2082     CContextClient* contextClientTmp ;
2083
2084     for (int i = 0; i < nbSrvPools; ++i)
2085     {
2086       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2087       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2088
2089       CEventClient event(getType(),EVENT_ID_PROCESS_GRID_ENABLED_FIELDS);
2090
2091       if (contextClientTmp->isServerLeader())
2092       {
2093         CMessage msg;
2094         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2095         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2096           event.push(*itRank,1,msg);
2097         contextClientTmp->sendEvent(event);
2098       }
2099       else contextClientTmp->sendEvent(event);
2100     }
2101   }
2102   CATCH_DUMP_ATTR
2103
2104   //! Server side: Receive a message to do some post processing
2105   void CContext::recvProcessingGridOfEnabledFields(CEventServer& event)
2106   TRY
2107   {
2108      CBufferIn* buffer=event.subEvents.begin()->buffer;
2109      // nothing to do, no call ??!!   
2110   }
2111   CATCH
2112
2113   //! Client side: Send a message to do some post processing on server
2114   void CContext::sendPostProcessing()
2115   TRY
2116   {
2117     int nbSrvPools ;
2118     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2119     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2120     else nbSrvPools = 0 ;
2121     CContextClient* contextClientTmp ;
2122
2123     for (int i = 0; i < nbSrvPools; ++i)
2124     {
2125       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2126       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2127       CEventClient event(getType(),EVENT_ID_POST_PROCESS);
2128       if (contextClientTmp->isServerLeader())
2129       {
2130         CMessage msg;
2131         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2132         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2133         event.push(*itRank,1,msg);
2134         contextClientTmp->sendEvent(event);
2135       }
2136       else contextClientTmp->sendEvent(event);
2137     }
2138   }
2139   CATCH_DUMP_ATTR
2140
2141   //! Server side: Receive a message to do some post processing
2142   void CContext::recvPostProcessing(CEventServer& event)
2143   TRY
2144   {
2145      CBufferIn* buffer=event.subEvents.begin()->buffer;
2146      getCurrent()->recvPostProcessing(*buffer);
2147   }
2148   CATCH
2149
2150   //! Server side: Receive a message to do some post processing
2151   void CContext::recvPostProcessing(CBufferIn& buffer)
2152   TRY
2153   {
2154      CCalendarWrapper::get(CCalendarWrapper::GetDefName())->createCalendar();
2155      postProcessing();
2156   }
2157   CATCH_DUMP_ATTR
2158
2159 
2160   /*!
2161   \brief Do some simple post processings after parsing xml file
2162      After the xml file (iodef.xml) is parsed, it is necessary to build all relations among
2163   created object, e.g: inhertance among fields, domain, axis. After that, all fiels as well as their parents (reference fields),
2164   which will be written out into netcdf files, are processed
2165   */
2166   void CContext::postProcessing()
2167   TRY
2168   {
2169     if (isPostProcessed) return;
2170
2171      // Make sure the calendar was correctly created
2172      if (!calendar)
2173        ERROR("CContext::postProcessing()", << "A calendar must be defined for the context \"" << getId() << "!\"")
2174      else if (calendar->getTimeStep() == NoneDu)
2175        ERROR("CContext::postProcessing()", << "A timestep must be defined for the context \"" << getId() << "!\"")
2176      // Calendar first update to set the current date equals to the start date
2177      calendar->update(0);
2178
2179      // Find all inheritance in xml structure
2180      this->solveAllInheritance();
2181
2182//      ShowTree(info(10));
2183
2184      // Check if some axis, domains or grids are eligible to for compressed indexed output.
2185      // Warning: This must be done after solving the inheritance and before the rest of post-processing
2186      checkAxisDomainsGridsEligibilityForCompressedOutput();      // only for field written on IO_SERVER service ????
2187
2188      // Check if some automatic time series should be generated
2189      // Warning: This must be done after solving the inheritance and before the rest of post-processing     
2190      prepareTimeseries();
2191
2192      //Initialisation du vecteur 'enabledFiles' contenant la liste des fichiers à sortir.
2193      findEnabledFiles();
2194      findEnabledWriteModeFiles();
2195      findEnabledReadModeFiles();
2196      findEnabledCouplerIn();
2197      findEnabledCouplerOut();
2198      createCouplerInterCommunicator() ;
2199
2200      // Find all enabled fields of each file     
2201      const vector<CField*>&& fileOutField = findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
2202      const vector<CField*>&& fileInField = findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
2203      const vector<CField*>&& CouplerOutField = findAllEnabledFieldsCouplerOut(this->enabledCouplerOut);
2204      const vector<CField*>&& CouplerInField = findAllEnabledFieldsCouplerIn(this->enabledCouplerIn);
2205
2206
2207
2208      // For now, only read files with client and only one level server
2209      // if (hasClient && !hasServer) findEnabledReadModeFiles();     
2210
2211
2212      // For now, only read files with client and only one level server
2213      // if (hasClient && !hasServer)
2214      //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
2215
2216      if (serviceType_==CServicesManager::CLIENT)
2217      {
2218        initReadFiles();
2219        // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
2220        this->readAttributesOfEnabledFieldsInReadModeFiles();
2221      }
2222
2223      // Only search and rebuild all reference objects of enable fields, don't transform
2224      this->solveOnlyRefOfEnabledFields();
2225
2226      // Search and rebuild all reference object of enabled fields, and transform
2227      this->solveAllRefOfEnabledFieldsAndTransform();
2228
2229      // Find all fields with read access from the public API
2230      if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
2231      // and solve the all reference for them
2232      if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
2233
2234      isPostProcessed = true;
2235   }
2236   CATCH_DUMP_ATTR
2237
2238   void CContext::createCouplerInterCommunicator(void)
2239   TRY
2240   {
2241      // juste for test now, in future need an scheduler to avoid dead-lock
2242      for(auto it=enabledCouplerOut.begin();it!=enabledCouplerOut.end();++it)
2243      {
2244        (*it)->createInterCommunicator() ;
2245      }
2246
2247      for(auto it=enabledCouplerIn.begin();it!=enabledCouplerIn.end();++it)
2248      {
2249        (*it)->createInterCommunicator() ;
2250      }
2251   }
2252   CATCH_DUMP_ATTR
2253
2254 
2255     //! Client side: Send infomation of active files (files are enabled to write out)
2256   void CContext::sendEnabledFiles(const std::vector<CFile*>& activeFiles)
2257   TRY
2258   {
2259     int size = activeFiles.size();
2260
2261     // In a context, each type has a root definition, e.g: axis, domain, field.
2262     // Every object must be a child of one of these root definition. In this case
2263     // all new file objects created on server must be children of the root "file_definition"
2264     StdString fileDefRoot("file_definition");
2265     CFileGroup* cfgrpPtr = CFileGroup::get(fileDefRoot);
2266
2267     for (int i = 0; i < size; ++i)
2268     {
2269       CFile* f = activeFiles[i];
2270       cfgrpPtr->sendCreateChild(f->getId(),f->getContextClient());
2271       f->sendAllAttributesToServer(f->getContextClient());
2272       f->sendAddAllVariables(f->getContextClient());
2273     }
2274   }
2275   CATCH_DUMP_ATTR
2276
2277   //! Client side: Send information of active fields (ones are written onto files)
2278   void CContext::sendEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
2279   TRY
2280   {
2281     int size = activeFiles.size();
2282     for (int i = 0; i < size; ++i)
2283     {
2284       activeFiles[i]->sendEnabledFields(activeFiles[i]->getContextClient());
2285     }
2286   }
2287   CATCH_DUMP_ATTR
2288
2289   //! Client side: Check if the defined axis, domains and grids are eligible for compressed indexed output
2290   void CContext::checkAxisDomainsGridsEligibilityForCompressedOutput()
2291   TRY
2292   {
2293     if (!(serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)) return;
2294
2295     const vector<CAxis*> allAxis = CAxis::getAll();
2296     for (vector<CAxis*>::const_iterator it = allAxis.begin(); it != allAxis.end(); it++)
2297       (*it)->checkEligibilityForCompressedOutput();
2298
2299     const vector<CDomain*> allDomains = CDomain::getAll();
2300     for (vector<CDomain*>::const_iterator it = allDomains.begin(); it != allDomains.end(); it++)
2301       (*it)->checkEligibilityForCompressedOutput();
2302
2303     const vector<CGrid*> allGrids = CGrid::getAll();
2304     for (vector<CGrid*>::const_iterator it = allGrids.begin(); it != allGrids.end(); it++)
2305       (*it)->checkEligibilityForCompressedOutput();
2306   }
2307   CATCH_DUMP_ATTR
2308
2309   //! Client side: Prepare the timeseries by adding the necessary files
2310   void CContext::prepareTimeseries()
2311   TRY
2312   {
2313     const std::vector<CFile*> allFiles = CFile::getAll();
2314     for (size_t i = 0; i < allFiles.size(); i++)
2315     {
2316       CFile* file = allFiles[i];
2317
2318       std::vector<CVariable*> fileVars, fieldVars, vars = file->getAllVariables();
2319       for (size_t k = 0; k < vars.size(); k++)
2320       {
2321         CVariable* var = vars[k];
2322
2323         if (var->ts_target.isEmpty()
2324              || var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both)
2325           fileVars.push_back(var);
2326
2327         if (!var->ts_target.isEmpty()
2328              && (var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both))
2329           fieldVars.push_back(var);
2330       }
2331
2332       if (!file->timeseries.isEmpty() && file->timeseries != CFile::timeseries_attr::none)
2333       {
2334         StdString fileNameStr("%file_name%") ;
2335         StdString tsPrefix = !file->ts_prefix.isEmpty() ? file->ts_prefix : fileNameStr ;
2336         
2337         StdString fileName=file->getFileOutputName();
2338         size_t pos=tsPrefix.find(fileNameStr) ;
2339         while (pos!=std::string::npos)
2340         {
2341           tsPrefix=tsPrefix.replace(pos,fileNameStr.size(),fileName) ;
2342           pos=tsPrefix.find(fileNameStr) ;
2343         }
2344       
2345         const std::vector<CField*> allFields = file->getAllFields();
2346         for (size_t j = 0; j < allFields.size(); j++)
2347         {
2348           CField* field = allFields[j];
2349
2350           if (!field->ts_enabled.isEmpty() && field->ts_enabled)
2351           {
2352             CFile* tsFile = CFile::create();
2353             tsFile->duplicateAttributes(file);
2354
2355             // Add variables originating from file and targeted to timeserie file
2356             for (size_t k = 0; k < fileVars.size(); k++)
2357               tsFile->getVirtualVariableGroup()->addChild(fileVars[k]);
2358
2359           
2360             tsFile->name = tsPrefix + "_";
2361             if (!field->name.isEmpty())
2362               tsFile->name.get() += field->name;
2363             else if (field->hasDirectFieldReference()) // We cannot use getBaseFieldReference() just yet
2364               tsFile->name.get() += field->field_ref;
2365             else
2366               tsFile->name.get() += field->getId();
2367
2368             if (!field->ts_split_freq.isEmpty())
2369               tsFile->split_freq = field->ts_split_freq;
2370
2371             CField* tsField = tsFile->addField();
2372             tsField->field_ref = field->getId();
2373
2374             // Add variables originating from file and targeted to timeserie field
2375             for (size_t k = 0; k < fieldVars.size(); k++)
2376               tsField->getVirtualVariableGroup()->addChild(fieldVars[k]);
2377
2378             vars = field->getAllVariables();
2379             for (size_t k = 0; k < vars.size(); k++)
2380             {
2381               CVariable* var = vars[k];
2382
2383               // Add variables originating from field and targeted to timeserie field
2384               if (var->ts_target.isEmpty()
2385                    || var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both)
2386                 tsField->getVirtualVariableGroup()->addChild(var);
2387
2388               // Add variables originating from field and targeted to timeserie file
2389               if (!var->ts_target.isEmpty()
2390                    && (var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both))
2391                 tsFile->getVirtualVariableGroup()->addChild(var);
2392             }
2393
2394             tsFile->solveFieldRefInheritance(true);
2395
2396             if (file->timeseries == CFile::timeseries_attr::exclusive)
2397               field->enabled = false;
2398           }
2399         }
2400
2401         // Finally disable the original file is need be
2402         if (file->timeseries == CFile::timeseries_attr::only)
2403          file->enabled = false;
2404       }
2405     }
2406   }
2407   CATCH_DUMP_ATTR
2408
2409   //! Client side: Send information of reference grid of active fields
2410   void CContext::sendRefGrid(const std::vector<CFile*>& activeFiles)
2411   TRY
2412   {
2413     std::set<pair<StdString,CContextClient*>> gridIds;
2414
2415     int sizeFile = activeFiles.size();
2416     CFile* filePtr(NULL);
2417
2418     // Firstly, find all reference grids of all active fields
2419     for (int i = 0; i < sizeFile; ++i)
2420     {
2421       filePtr = activeFiles[i];
2422       std::vector<CField*> enabledFields = filePtr->getEnabledFields();
2423       int sizeField = enabledFields.size();
2424       for (int numField = 0; numField < sizeField; ++numField)
2425       {
2426         if (0 != enabledFields[numField]->getRelGrid())
2427           gridIds.insert(make_pair(CGrid::get(enabledFields[numField]->getRelGrid())->getId(),enabledFields[numField]->getContextClient()));
2428       }
2429     }
2430
2431     // Create all reference grids on server side
2432     StdString gridDefRoot("grid_definition");
2433     CGridGroup* gridPtr = CGridGroup::get(gridDefRoot);
2434     for (auto it = gridIds.begin(); it != gridIds.end(); ++it)
2435     {
2436       gridPtr->sendCreateChild(it->first,it->second);
2437       CGrid::get(it->first)->sendAllAttributesToServer(it->second);
2438       CGrid::get(it->first)->sendAllDomains(it->second);
2439       CGrid::get(it->first)->sendAllAxis(it->second);
2440       CGrid::get(it->first)->sendAllScalars(it->second);
2441     }
2442   }
2443   CATCH_DUMP_ATTR
2444
2445   //! Client side: Send information of reference domain, axis and scalar of active fields
2446   void CContext::sendRefDomainsAxisScalars(const std::vector<CFile*>& activeFiles)
2447   TRY
2448   {
2449     std::set<pair<StdString,CContextClient*>> domainIds, axisIds, scalarIds;
2450
2451     // Find all reference domain and axis of all active fields
2452     int numEnabledFiles = activeFiles.size();
2453     for (int i = 0; i < numEnabledFiles; ++i)
2454     {
2455       std::vector<CField*> enabledFields = activeFiles[i]->getEnabledFields();
2456       int numEnabledFields = enabledFields.size();
2457       for (int j = 0; j < numEnabledFields; ++j)
2458       {
2459         CContextClient* contextClient=enabledFields[j]->getContextClient() ;
2460         const std::vector<StdString>& prDomAxisScalarId = enabledFields[j]->getRefDomainAxisIds();
2461         if ("" != prDomAxisScalarId[0]) domainIds.insert(make_pair(prDomAxisScalarId[0],contextClient));
2462         if ("" != prDomAxisScalarId[1]) axisIds.insert(make_pair(prDomAxisScalarId[1],contextClient));
2463         if ("" != prDomAxisScalarId[2]) scalarIds.insert(make_pair(prDomAxisScalarId[2],contextClient));
2464       }
2465     }
2466
2467     // Create all reference axis on server side
2468     std::set<StdString>::iterator itDom, itAxis, itScalar;
2469     std::set<StdString>::const_iterator itE;
2470
2471     StdString scalarDefRoot("scalar_definition");
2472     CScalarGroup* scalarPtr = CScalarGroup::get(scalarDefRoot);
2473     
2474     for (auto itScalar = scalarIds.begin(); itScalar != scalarIds.end(); ++itScalar)
2475     {
2476       if (!itScalar->first.empty())
2477       {
2478         scalarPtr->sendCreateChild(itScalar->first,itScalar->second);
2479         CScalar::get(itScalar->first)->sendAllAttributesToServer(itScalar->second);
2480       }
2481     }
2482
2483     StdString axiDefRoot("axis_definition");
2484     CAxisGroup* axisPtr = CAxisGroup::get(axiDefRoot);
2485     
2486     for (auto itAxis = axisIds.begin(); itAxis != axisIds.end(); ++itAxis)
2487     {
2488       if (!itAxis->first.empty())
2489       {
2490         axisPtr->sendCreateChild(itAxis->first, itAxis->second);
2491         CAxis::get(itAxis->first)->sendAllAttributesToServer(itAxis->second);
2492       }
2493     }
2494
2495     // Create all reference domains on server side
2496     StdString domDefRoot("domain_definition");
2497     CDomainGroup* domPtr = CDomainGroup::get(domDefRoot);
2498     
2499     for (auto itDom = domainIds.begin(); itDom != domainIds.end(); ++itDom)
2500     {
2501       if (!itDom->first.empty()) {
2502          domPtr->sendCreateChild(itDom->first, itDom->second);
2503          CDomain::get(itDom->first)->sendAllAttributesToServer(itDom->second);
2504       }
2505     }
2506   }
2507   CATCH_DUMP_ATTR
2508
2509   //! Update calendar in each time step
2510   void CContext::updateCalendar(int step)
2511   TRY
2512   {
2513      int prevStep = calendar->getStep();
2514
2515      if (prevStep < step)
2516      {
2517        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2518        {
2519          doPreTimestepOperationsForEnabledReadModeFiles();
2520        }
2521
2522        info(50) << "updateCalendar : before : " << calendar->getCurrentDate() << endl;
2523        calendar->update(step);
2524        info(50) << "updateCalendar : after : " << calendar->getCurrentDate() << endl;
2525  #ifdef XIOS_MEMTRACK_LIGHT
2526        info(50) << " Current memory used by XIOS : "<<  MemTrack::getCurrentMemorySize()*1.0/(1024*1024)<<" Mbyte, at timestep "<<step<<" of context "<<this->getId()<<endl ;
2527  #endif
2528
2529        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2530        {
2531          doPostTimestepOperationsForEnabledReadModeFiles();
2532          garbageCollector.invalidate(calendar->getCurrentDate());
2533        }
2534      }
2535      else if (prevStep == step)
2536        info(50) << "updateCalendar: already at step " << step << ", no operation done." << endl;
2537      else // if (prevStep > step)
2538        ERROR("void CContext::updateCalendar(int step)",
2539              << "Illegal calendar update: previous step was " << prevStep << ", new step " << step << "is in the past!")
2540   }
2541   CATCH_DUMP_ATTR
2542
2543   void CContext::initReadFiles(void)
2544   TRY
2545   {
2546      vector<CFile*>::const_iterator it;
2547
2548      for (it=enabledReadModeFiles.begin(); it != enabledReadModeFiles.end(); it++)
2549      {
2550         (*it)->initRead();
2551      }
2552   }
2553   CATCH_DUMP_ATTR
2554
2555   //! Server side: Create header of netcdf file
2556   void CContext::createFileHeader(void)
2557   TRY
2558   {
2559      vector<CFile*>::const_iterator it;
2560
2561      for (it=enabledFiles.begin(); it != enabledFiles.end(); it++)
2562      // for (it=enabledWriteModeFiles.begin(); it != enabledWriteModeFiles.end(); it++)
2563      {
2564         (*it)->initWrite();
2565      }
2566   }
2567   CATCH_DUMP_ATTR
2568
2569   //! Get current context
2570   CContext* CContext::getCurrent(void)
2571   TRY
2572   {
2573     return CObjectFactory::GetObject<CContext>(CObjectFactory::GetCurrentContextId()).get();
2574   }
2575   CATCH
2576
2577   /*!
2578   \brief Set context with an id be the current context
2579   \param [in] id identity of context to be set to current
2580   */
2581   void CContext::setCurrent(const string& id)
2582   TRY
2583   {
2584     CObjectFactory::SetCurrentContextId(id);
2585     CGroupFactory::SetCurrentContextId(id);
2586   }
2587   CATCH
2588
2589  /*!
2590  \brief Create a context with specific id
2591  \param [in] id identity of new context
2592  \return pointer to the new context or already-existed one with identity id
2593  */
2594  CContext* CContext::create(const StdString& id)
2595  TRY
2596  {
2597    CContext::setCurrent(id);
2598
2599    bool hasctxt = CContext::has(id);
2600    CContext* context = CObjectFactory::CreateObject<CContext>(id).get();
2601    getRoot();
2602    if (!hasctxt) CGroupFactory::AddChild(root, context->getShared());
2603
2604#define DECLARE_NODE(Name_, name_) \
2605    C##Name_##Definition::create(C##Name_##Definition::GetDefName());
2606#define DECLARE_NODE_PAR(Name_, name_)
2607#include "node_type.conf"
2608
2609    return (context);
2610  }
2611  CATCH
2612
2613     //! Server side: Receive a message to do some post processing
2614  void CContext::recvRegistry(CEventServer& event)
2615  TRY
2616  {
2617    CBufferIn* buffer=event.subEvents.begin()->buffer;
2618    getCurrent()->recvRegistry(*buffer);
2619  }
2620  CATCH
2621
2622  void CContext::recvRegistry(CBufferIn& buffer)
2623  TRY
2624  {
2625    if (server->intraCommRank==0)
2626    {
2627      CRegistry registry(server->intraComm) ;
2628      registry.fromBuffer(buffer) ;
2629      registryOut->mergeRegistry(registry) ;
2630    }
2631  }
2632  CATCH_DUMP_ATTR
2633
2634  void CContext::sendRegistry(void)
2635  TRY
2636  {
2637    registryOut->hierarchicalGatherRegistry() ;
2638
2639    int nbSrvPools ;
2640    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2641    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2642    else nbSrvPools = 0 ;
2643    CContextClient* contextClientTmp ;
2644
2645    for (int i = 0; i < nbSrvPools; ++i)
2646    {
2647      if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2648      else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2649
2650      CEventClient event(CContext::GetType(), CContext::EVENT_ID_SEND_REGISTRY);
2651      if (contextClientTmp->isServerLeader())
2652      {
2653        CMessage msg ;
2654        if (contextClientTmp->clientRank==0) msg<<*registryOut ;
2655        const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2656        for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2657             event.push(*itRank,1,msg);
2658        contextClientTmp->sendEvent(event);
2659      }
2660      else contextClientTmp->sendEvent(event);
2661    }
2662  }
2663  CATCH_DUMP_ATTR
2664
2665 
2666  void CContext::sendFinalizeClient(CContextClient* contextClient, const string& contextClientId)
2667  TRY
2668  {
2669    CEventClient event(getType(),EVENT_ID_CONTEXT_FINALIZE_CLIENT);
2670    if (contextClient->isServerLeader())
2671    {
2672      CMessage msg;
2673      msg<<contextClientId ;
2674      const std::list<int>& ranks = contextClient->getRanksServerLeader();
2675      for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2676           event.push(*itRank,1,msg);
2677      contextClient->sendEvent(event);
2678    }
2679    else contextClient->sendEvent(event);
2680  }
2681  CATCH_DUMP_ATTR
2682
2683 
2684  void CContext::recvFinalizeClient(CEventServer& event)
2685  TRY
2686  {
2687    CBufferIn* buffer=event.subEvents.begin()->buffer;
2688    string id;
2689    *buffer>>id;
2690    get(id)->recvFinalizeClient(*buffer);
2691  }
2692  CATCH
2693
2694  void CContext::recvFinalizeClient(CBufferIn& buffer)
2695  TRY
2696  {
2697    countChildContextFinalized_++ ;
2698  }
2699  CATCH_DUMP_ATTR
2700
2701
2702
2703
2704  /*!
2705  * \fn bool CContext::isFinalized(void)
2706  * Context is finalized if it received context post finalize event.
2707  */
2708  bool CContext::isFinalized(void)
2709  TRY
2710  {
2711    return finalized;
2712  }
2713  CATCH_DUMP_ATTR
2714  ///--------------------------------------------------------------
2715  StdString CContext::dumpClassAttributes(void)
2716  {
2717    StdString str;
2718    str.append("enabled files=\"");
2719    int size = this->enabledFiles.size();
2720    for (int i = 0; i < size; ++i)
2721    {
2722      str.append(enabledFiles[i]->getId());
2723      str.append(" ");
2724    }
2725    str.append("\"");
2726    return str;
2727  }
2728
2729} // namespace xios
Note: See TracBrowser for help on using the repository browser.