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

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

Big update on on going work related to data distribution and transfer between clients and servers.
Revisite of the source and store filter using "connectors".

-> inputs work again

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: 109.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 couplerOut : couplerOutClient_)
711      if (!finalized) couplerOut.second->checkBuffers();
712   
713    for (auto couplerIn : couplerInClient_)
714      if (!finalized) couplerIn.second->checkBuffers();
715   
716    for (auto couplerOut : couplerOutServer_)
717      if (!finalized) couplerOut.second->eventLoop(enableEventsProcessing);
718
719    for (auto couplerIn : couplerInServer_)
720      if (!finalized) couplerIn.second->eventLoop(enableEventsProcessing);
721   
722    if (server!=nullptr) if (!finalized) finished &= server->eventLoop(enableEventsProcessing);
723 
724    return finalized && finished ;
725  }
726
727  void CContext::addCouplingChanel(const std::string& fullContextId, bool out)
728  {
729     int contextLeader ;
730     
731     if (out)
732     { 
733       if (couplerOutClient_.find(fullContextId)==couplerOutClient_.end()) 
734       {
735         bool ok=CXios::getContextsManager()->getContextLeader(fullContextId, contextLeader, getIntraComm()) ;
736     
737         MPI_Comm interComm, interCommClient, interCommServer  ;
738         MPI_Comm intraCommClient, intraCommServer ;
739
740         if (ok) MPI_Intercomm_create(getIntraComm(), 0, CXios::getXiosComm(), contextLeader, 0, &interComm) ;
741
742        MPI_Comm_dup(intraComm_, &intraCommClient) ;
743        MPI_Comm_dup(intraComm_, &intraCommServer) ;
744        MPI_Comm_dup(interComm, &interCommClient) ;
745        MPI_Comm_dup(interComm, &interCommServer) ;
746        CContextClient* client = new CContextClient(this, intraCommClient, interCommClient);
747        CContextServer* server = new CContextServer(this, intraCommServer, interCommServer);
748        client->setAssociatedServer(server) ;
749        server->setAssociatedClient(client) ;
750        MPI_Comm_free(&interComm) ;
751        couplerOutClient_[fullContextId] = client ;
752        couplerOutServer_[fullContextId] = server ;
753
754/*
755      // for now, we don't now which beffer size must be used for client coupler
756      // It will be evaluated later. Fix a constant size for now...
757      // set to 10Mb for development
758       map<int,size_t> bufferSize, maxEventSize ;
759       for(int i=0;i<client->getRemoteSize();i++)
760       {
761         bufferSize[i]=10000000 ;
762         maxEventSize[i]=10000000 ;
763       }
764
765       client->setBufferSize(bufferSize, maxEventSize);   
766*/
767      }
768    }
769    else if (couplerInClient_.find(fullContextId)==couplerInClient_.end())
770    {
771      bool ok=CXios::getContextsManager()->getContextLeader(fullContextId, contextLeader, getIntraComm()) ;
772     
773       MPI_Comm interComm, interCommClient, interCommServer  ;
774       MPI_Comm intraCommClient, intraCommServer ;
775
776       if (ok) MPI_Intercomm_create(getIntraComm(), 0, CXios::getXiosComm(), contextLeader, 0, &interComm) ;
777
778       MPI_Comm_dup(intraComm_, &intraCommClient) ;
779       MPI_Comm_dup(intraComm_, &intraCommServer) ;
780       MPI_Comm_dup(interComm, &interCommServer) ;
781       MPI_Comm_dup(interComm, &interCommClient) ;
782       CContextServer* server = new CContextServer(this, intraCommServer, interCommServer);
783       CContextClient* client = new CContextClient(this, intraCommClient, interCommClient);
784       client->setAssociatedServer(server) ;
785       server->setAssociatedClient(client) ;
786       MPI_Comm_free(&interComm) ;
787
788       map<int,size_t> bufferSize, maxEventSize ;
789       for(int i=0;i<client->getRemoteSize();i++)
790       {
791         bufferSize[i]=10000000 ;
792         maxEventSize[i]=10000000 ;
793       }
794
795       client->setBufferSize(bufferSize, maxEventSize);   
796       couplerInClient_[fullContextId] = client ;
797       couplerInServer_[fullContextId] = server ;       
798    }
799  }
800 
801  void CContext::globalEventLoop(void)
802  {
803    CXios::getDaemonsManager()->eventLoop() ;
804    setCurrent(getId()) ;
805  }
806
807
808   void CContext::finalize(void)
809   TRY
810   {
811      registryOut->hierarchicalGatherRegistry() ;
812      if (server->intraCommRank==0) CXios::globalRegistry->mergeRegistry(*registryOut) ;
813
814      if (serviceType_==CServicesManager::CLIENT)
815      {
816//ym        doPreTimestepOperationsForEnabledReadModeFiles(); // For now we only use server level 1 to read data
817
818        triggerLateFields() ;
819
820        // inform couplerIn that I am finished
821        for(auto& couplerInClient : couplerInClient_) sendCouplerInContextFinalized(couplerInClient.second) ;
822
823        // wait until received message from couplerOut that they have finished
824        bool couplersInFinalized ;
825        do
826        {
827          couplersInFinalized=true ;
828          for(auto& couplerOutClient : couplerOutClient_) couplersInFinalized &= isCouplerInContextFinalized(couplerOutClient.second) ; 
829          globalEventLoop() ;
830        } while (!couplersInFinalized) ;
831
832        info(100)<<"DEBUG: context "<<getId()<<" Send client finalize"<<endl ;
833        client->finalize();
834        info(100)<<"DEBUG: context "<<getId()<<" Client finalize sent"<<endl ;
835        while (client->havePendingRequests()) client->checkBuffers();
836        info(100)<<"DEBUG: context "<<getId()<<" no pending request ok"<<endl ;
837        bool notifiedFinalized=false ;
838        do
839        {
840          notifiedFinalized=client->isNotifiedFinalized() ;
841        } while (!notifiedFinalized) ;
842        client->releaseBuffers();
843        info(100)<<"DEBUG: context "<<getId()<<" release client ok"<<endl ;
844      }
845      else if (serviceType_==CServicesManager::GATHERER)
846      {
847         for (int i = 0; i < clientPrimServer.size(); ++i)
848         {
849           clientPrimServer[i]->finalize();
850           bool bufferReleased;
851           do
852           {
853             clientPrimServer[i]->checkBuffers();
854             bufferReleased = !clientPrimServer[i]->havePendingRequests();
855           } while (!bufferReleased);
856           
857           bool notifiedFinalized=false ;
858           do
859           {
860             notifiedFinalized=clientPrimServer[i]->isNotifiedFinalized() ;
861           } while (!notifiedFinalized) ;
862           clientPrimServer[i]->releaseBuffers();
863         }
864         closeAllFile();
865
866      }
867      else if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)
868      {
869        closeAllFile();
870      }
871
872      freeComms() ;
873       
874      parentServerContext_->freeComm() ;
875      finalized = true;
876      info(20)<<"CContext: Context <"<<getId()<<"> is finalized."<<endl;
877   }
878   CATCH_DUMP_ATTR
879
880   //! Free internally allocated communicators
881   void CContext::freeComms(void)
882   TRY
883   {
884     for (std::list<MPI_Comm>::iterator it = comms.begin(); it != comms.end(); ++it)
885       MPI_Comm_free(&(*it));
886     comms.clear();
887   }
888   CATCH_DUMP_ATTR
889
890   //! Deallocate buffers allocated by clientContexts
891   void CContext::releaseClientBuffers(void)
892   TRY
893   {
894     client->releaseBuffers();
895     for (int i = 0; i < clientPrimServer.size(); ++i)
896       clientPrimServer[i]->releaseBuffers();
897   }
898   CATCH_DUMP_ATTR
899
900   void CContext::postProcessingGlobalAttributes()
901   TRY
902   {
903     if (allProcessed) return; 
904     
905    // create intercommunicator with servers.
906    // not sure it is the good place to be called here
907    createServerInterComm() ;
908
909
910     // After xml is parsed, there are some more works with post processing
911     postProcessing();
912
913     // Distribute files between secondary servers according to the data size
914     distributeFiles(this->enabledWriteModeFiles);
915
916     // Check grid and calculate its distribution
917     checkGridEnabledFields();
918
919     setClientServerBuffer(client, (serviceType_==CServicesManager::CLIENT) ) ;
920     for (int i = 0; i < clientPrimServer.size(); ++i)
921         setClientServerBuffer(clientPrimServer[i], true);
922
923   
924     if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
925     { 
926       if (serviceType_==CServicesManager::GATHERER)
927       { 
928         for (auto it=clientPrimServer.begin(); it!=clientPrimServer.end();++it) 
929         {
930           this->sendAllAttributesToServer(*it); // Send all attributes of current context to server
931           CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(*it); // Send all attributes of current calendar
932         }
933       }
934       else 
935       {
936         this->sendAllAttributesToServer(client);   // Send all attributes of current context to server
937         CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current calendar
938       }
939
940
941      // We have enough information to send to server
942      // First of all, send all enabled files
943      sendEnabledFiles(this->enabledWriteModeFiles);
944      // We only use server-level 1 (for now) to read data
945      if (serviceType_==CServicesManager::CLIENT)  sendEnabledFiles(this->enabledReadModeFiles);
946
947      // Then, send all enabled fields     
948      sendEnabledFieldsInFiles(this->enabledWriteModeFiles);
949     
950      if (serviceType_==CServicesManager::CLIENT) sendEnabledFieldsInFiles(this->enabledReadModeFiles);
951
952      // Then, check whether we have domain_ref, axis_ref or scalar_ref attached to the enabled fields
953      // If any, so send them to server
954       sendRefDomainsAxisScalars(this->enabledWriteModeFiles); 
955     
956      if (serviceType_==CServicesManager::CLIENT) sendRefDomainsAxisScalars(this->enabledReadModeFiles);       
957
958       // Check whether enabled fields have grid_ref, if any, send this info to server
959      sendRefGrid(this->enabledFiles);
960      // This code may be useful in the future when we want to seperate completely read and write
961      // sendRefGrid(this->enabledWriteModeFiles);
962      // if (!hasServer)
963      //   sendRefGrid(this->enabledReadModeFiles);
964     
965      // A grid of enabled fields composed of several components which must be checked then their
966      // checked attributes should be sent to server
967      sendGridComponentEnabledFieldsInFiles(this->enabledFiles); // This code can be seperated in two (one for reading, another for writing)
968
969       // We have a xml tree on the server side and now, it should be also processed
970      sendPostProcessing();
971       
972      // Finally, we send information of grid itself to server
973      sendGridEnabledFieldsInFiles(this->enabledWriteModeFiles);       
974     
975      if (serviceType_==CServicesManager::CLIENT) sendGridEnabledFieldsInFiles(this->enabledReadModeFiles);       
976     
977     }
978     allProcessed = true;
979   }
980   CATCH_DUMP_ATTR
981
982   void CContext::sendPostProcessingGlobalAttributes()
983   TRY
984   {
985
986    int nbSrvPools ;
987    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
988    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
989    else nbSrvPools = 0 ;
990    CContextClient* contextClientTmp ;
991
992    for (int i = 0; i < nbSrvPools; ++i)
993     {
994       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
995       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
996     
997       CEventClient event(getType(),EVENT_ID_POST_PROCESS_GLOBAL_ATTRIBUTES);
998
999       if (contextClientTmp->isServerLeader())
1000       {
1001         CMessage msg;
1002         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
1003         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1004           event.push(*itRank,1,msg);
1005         contextClientTmp->sendEvent(event);
1006       }
1007       else contextClientTmp->sendEvent(event);
1008     }
1009   }
1010   CATCH_DUMP_ATTR
1011
1012   void CContext::recvPostProcessingGlobalAttributes(CEventServer& event)
1013   TRY
1014   {
1015      CBufferIn* buffer=event.subEvents.begin()->buffer;
1016      getCurrent()->recvPostProcessingGlobalAttributes(*buffer);
1017   }
1018   CATCH
1019
1020   void CContext::recvPostProcessingGlobalAttributes(CBufferIn& buffer)
1021   TRY
1022   {     
1023      postProcessingGlobalAttributes();
1024   }
1025   CATCH_DUMP_ATTR
1026
1027   /*!
1028   \brief Close all the context defintion and do processing data
1029      After everything is well defined on client side, they will be processed and sent to server
1030   From the version 2.0, sever and client work no more on the same database. Moreover, client(s) will send
1031   all necessary information to server, from which each server can build its own database.
1032   Because the role of server is to write out field data on a specific netcdf file,
1033   the only information that it needs is the enabled files
1034   and the active fields (fields will be written onto active files)
1035   */
1036  void CContext::closeDefinition(void)
1037   TRY
1038   {
1039     CTimer::get("Context : close definition").resume() ;
1040     
1041     // create intercommunicator with servers.
1042     // not sure it is the good place to be called here
1043     createServerInterComm() ;
1044
1045
1046     // After xml is parsed, there are some more works with post processing
1047//     postProcessing();
1048
1049   
1050    // Make sure the calendar was correctly created
1051    if (serviceType_!=CServicesManager::CLIENT) CCalendarWrapper::get(CCalendarWrapper::GetDefName())->createCalendar();
1052    if (!calendar)
1053      ERROR("CContext::postProcessing()", << "A calendar must be defined for the context \"" << getId() << "!\"")
1054    else if (calendar->getTimeStep() == NoneDu)
1055      ERROR("CContext::postProcessing()", << "A timestep must be defined for the context \"" << getId() << "!\"")
1056    // Calendar first update to set the current date equals to the start date
1057    calendar->update(0);
1058
1059    // Résolution des héritages descendants (càd des héritages de groupes)
1060    // pour chacun des contextes.
1061    solveDescInheritance(true);
1062 
1063    // Solve inheritance for field to know if enabled or not.
1064    for (auto field : CField::getAll()) field->solveRefInheritance();
1065
1066    // Check if some axis, domains or grids are eligible to for compressed indexed output.
1067    // Warning: This must be done after solving the inheritance and before the rest of post-processing
1068    // --> later ????    checkAxisDomainsGridsEligibilityForCompressedOutput();     
1069
1070      // Check if some automatic time series should be generated
1071      // Warning: This must be done after solving the inheritance and before the rest of post-processing     
1072
1073    // The timeseries should only be prepared in client
1074    prepareTimeseries();
1075
1076    //Initialisation du vecteur 'enabledFiles' contenant la liste des fichiers à sortir.
1077    findEnabledFiles();
1078    findEnabledWriteModeFiles();
1079    findEnabledReadModeFiles();
1080    findEnabledCouplerIn();
1081    findEnabledCouplerOut();
1082    createCouplerInterCommunicator() ;
1083
1084    // Find all enabled fields of each file     
1085    vector<CField*>&& fileOutField = findAllEnabledFieldsInFileOut(this->enabledWriteModeFiles);
1086    vector<CField*>&& fileInField = findAllEnabledFieldsInFileIn(this->enabledReadModeFiles);
1087    vector<CField*>&& couplerOutField = findAllEnabledFieldsCouplerOut(this->enabledCouplerOut);
1088    vector<CField*>&& couplerInField = findAllEnabledFieldsCouplerIn(this->enabledCouplerIn);
1089    findFieldsWithReadAccess();
1090    vector<CField*>& fieldWithReadAccess = fieldsWithReadAccess_ ;
1091    vector<CField*> fieldModelIn ; // fields potentially from model
1092     
1093    // define if files are on clientSied or serverSide
1094    if (serviceType_==CServicesManager::CLIENT)
1095    {
1096      for (auto& file : enabledWriteModeFiles) file->setClientSide() ;
1097      for (auto& file : enabledReadModeFiles) file->setClientSide() ;
1098    }
1099    else
1100    {
1101      for (auto& file : enabledWriteModeFiles) file->setServerSide() ;
1102      for (auto& file : enabledReadModeFiles) file->setServerSide() ;
1103    }
1104
1105   
1106    for (auto& field : couplerInField)
1107    {
1108      field->unsetGridCompleted() ;
1109    }
1110// find all field potentially at workflow end
1111    vector<CField*> endWorkflowFields ;
1112    endWorkflowFields.reserve(fileOutField.size()+couplerOutField.size()+fieldWithReadAccess.size()) ;
1113    endWorkflowFields.insert(endWorkflowFields.end(),fileOutField.begin(), fileOutField.end()) ;
1114    endWorkflowFields.insert(endWorkflowFields.end(),couplerOutField.begin(), couplerOutField.end()) ;
1115    endWorkflowFields.insert(endWorkflowFields.end(),fieldWithReadAccess.begin(), fieldWithReadAccess.end()) ;
1116
1117    bool workflowGraphIsCompleted ;
1118
1119    bool first=true ;
1120    do
1121    {
1122      workflowGraphIsCompleted=true; 
1123      for(auto endWorkflowField : endWorkflowFields) 
1124      {
1125        workflowGraphIsCompleted &= endWorkflowField->buildWorkflowGraph(garbageCollector) ;
1126      }
1127
1128      for(auto couplerIn : enabledCouplerIn) couplerIn->assignContext() ;
1129      for(auto field : couplerInField) field->makeGridAliasForCoupling();
1130      for(auto field : couplerInField) this->sendCouplerInReady(field->getContextClient()) ;
1131   
1132
1133      // assign context to coupler out and related fields
1134      for(auto couplerOut : enabledCouplerOut) couplerOut->assignContext() ;
1135      // for now supose that all coupling out endpoint are succesfull. The difficultie is client/server buffer evaluation
1136      for(auto field : couplerOutField) 
1137      {
1138        field->computeGridIndexToFileServer() ; // same kind of index than for file server -> in future distribution may change
1139      }
1140      if (first) setClientServerBuffer(couplerOutField, true) ; // set buffer context
1141
1142      bool couplersReady ;
1143      do 
1144      {
1145        couplersReady=true ;
1146        for(auto field : couplerOutField)
1147        {
1148          bool ready = isCouplerInReady(field->getContextClient()) ; 
1149          if (ready) field->sendFieldToCouplerOut() ;
1150          couplersReady &= ready ;
1151        }
1152        if (!couplersReady) this->eventLoop() ;
1153      } while (!couplersReady) ;
1154
1155      first=false ;
1156      this->eventLoop() ;
1157    } while (!workflowGraphIsCompleted) ;
1158
1159    for( auto field : couplerInField) couplerInFields_.push_back(field) ;
1160
1161    // get all field coming potentially from model
1162    for (auto field : CField::getAll() ) if (field->getModelIn()) fieldModelIn.push_back(field) ;
1163
1164    // Distribute files between secondary servers according to the data size => assign a context to a file and then to fields
1165    if (serviceType_==CServicesManager::GATHERER) distributeFiles(this->enabledWriteModeFiles);
1166    else if (serviceType_==CServicesManager::CLIENT) for(auto file : this->enabledWriteModeFiles) file->setContextClient(client) ;
1167
1168    // client side, assign context for file reading
1169    if (serviceType_==CServicesManager::CLIENT) for(auto file : this->enabledReadModeFiles) file->setContextClient(client) ;
1170   
1171    // server side, assign context where to send file data read
1172    if (serviceType_==CServicesManager::CServicesManager::GATHERER || serviceType_==CServicesManager::IO_SERVER) 
1173      for(auto file : this->enabledReadModeFiles) file->setContextClient(client) ;
1174   
1175    // workflow endpoint => sent to IO/SERVER
1176    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1177    {
1178      for(auto field : fileOutField) 
1179      {
1180        field->connectToFileServer(garbageCollector) ; // connect the field to server filter
1181        field->computeGridIndexToFileServer() ; // compute grid index for transfer to the server context
1182      }
1183      setClientServerBuffer(fileOutField, true) ; // set buffer context
1184      for(auto field : fileOutField) field->sendFieldToFileServer() ;
1185    }
1186
1187    // workflow endpoint => write to file
1188    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)
1189    {
1190      for(auto field : fileOutField) 
1191      {
1192        field->connectToFileWriter(garbageCollector) ; // connect the field to server filter
1193      }
1194    }
1195   
1196    // workflow endpoint => Send data from server to client
1197    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER)
1198    {
1199      for(auto field : fileInField) 
1200      {
1201        field->connectToServerToClient(garbageCollector) ;
1202      }
1203    }
1204
1205    // workflow endpoint => sent to model on client side
1206    if (serviceType_==CServicesManager::CLIENT)
1207    {
1208      for(auto field : fieldWithReadAccess) field->connectToModelOutput(garbageCollector) ;
1209    }
1210
1211
1212    // workflow startpoint => data from model
1213    if (serviceType_==CServicesManager::CLIENT)
1214    {
1215      for(auto field : fieldModelIn) 
1216      {
1217        field->connectToModelInput(garbageCollector) ; // connect the field to server filter
1218        // grid index will be computed on the fly
1219      }
1220    }
1221   
1222    // workflow startpoint => data from client on server side
1223    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER || serviceType_==CServicesManager::OUT_SERVER)
1224    {
1225      for(auto field : fieldModelIn) 
1226      {
1227        field->connectToClientInput(garbageCollector) ; // connect the field to server filter
1228      }
1229    }
1230
1231   
1232    for(auto field : couplerInField) 
1233    {
1234      field->connectToCouplerIn(garbageCollector) ; // connect the field to server filter
1235    }
1236   
1237   
1238    for(auto field : couplerOutField) 
1239    {
1240      field->connectToCouplerOut(garbageCollector) ; // for now the same kind of filter that for file server
1241    }
1242
1243     // workflow startpoint => data from server on client side
1244    if (serviceType_==CServicesManager::CLIENT)
1245    {
1246      for(auto field : fileInField) 
1247      {
1248        field->connectToServerInput(garbageCollector) ; // connect the field to server filter
1249    // obsolete    field->computeGridIndexToFileServer() ; // compute grid index for transfer to the server context
1250        field->sendFieldToInputFileServer() ;
1251        fileInFields_.push_back(field) ;
1252      }
1253    }
1254
1255    // workflow startpoint => data read from file on server side
1256    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER)
1257    {
1258      for(auto field : fileInField) 
1259      {
1260        field->connectToFileReader(garbageCollector) ;
1261      }
1262    }
1263   
1264    // construct slave server list
1265    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) 
1266    {
1267      for(auto field : fileOutField) slaveServers_.insert(field->getContextClient()) ; 
1268      for(auto field : fileInField)  slaveServers_.insert(field->getContextClient()) ; 
1269    }
1270
1271    for(auto& slaveServer : slaveServers_) sendCloseDefinition(slaveServer) ;
1272
1273    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER) 
1274    {
1275      createFileHeader();
1276    }
1277
1278    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1279   
1280    // send signal to couplerIn context that definition phasis is done
1281
1282    for(auto& couplerInClient : couplerInClient_) sendCouplerInCloseDefinition(couplerInClient.second) ;
1283
1284    // wait until all couplerIn signal that closeDefition is done.
1285    bool ok;
1286    do
1287    {
1288      ok = true ;
1289      for(auto& couplerOutClient : couplerOutClient_) ok &= isCouplerInCloseDefinition(couplerOutClient.second) ;
1290      this->eventLoop() ; 
1291    } while (!ok) ;
1292
1293    return ;
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304    // For now, only read files with client and only one level server
1305    // if (hasClient && !hasServer) findEnabledReadModeFiles();     
1306
1307    // Find all enabled fields of each file     
1308    findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
1309    findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
1310
1311    // For now, only read files with client and only one level server
1312    // if (hasClient && !hasServer)
1313    //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
1314
1315    if (serviceType_==CServicesManager::CLIENT)
1316    {
1317      initReadFiles();
1318      // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
1319      this->readAttributesOfEnabledFieldsInReadModeFiles();
1320    }
1321
1322    // Only search and rebuild all reference objects of enable fields, don't transform
1323    this->solveOnlyRefOfEnabledFields();
1324
1325    // Search and rebuild all reference object of enabled fields, and transform
1326    this->solveAllRefOfEnabledFieldsAndTransform();
1327
1328    // Find all fields with read access from the public API
1329    if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
1330    // and solve the all reference for them
1331    if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
1332
1333    isPostProcessed = true;
1334
1335
1336
1337    // Distribute files between secondary servers according to the data size
1338    distributeFiles(this->enabledWriteModeFiles);
1339
1340    // Check grid and calculate its distribution
1341    checkGridEnabledFields();
1342
1343    setClientServerBuffer(client, (serviceType_==CServicesManager::CLIENT) ) ;
1344    for (int i = 0; i < clientPrimServer.size(); ++i)
1345         setClientServerBuffer(clientPrimServer[i], true);
1346
1347   
1348    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1349    { 
1350      if (serviceType_==CServicesManager::GATHERER)
1351      { 
1352        for (auto it=clientPrimServer.begin(); it!=clientPrimServer.end();++it) 
1353        {
1354          this->sendAllAttributesToServer(*it); // Send all attributes of current context to server
1355          CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(*it); // Send all attributes of current calendar
1356        }
1357      }
1358      else 
1359      {
1360        this->sendAllAttributesToServer(client);   // Send all attributes of current context to server
1361        CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current calendar
1362      }
1363
1364      // We have enough information to send to server
1365      // First of all, send all enabled files
1366      sendEnabledFiles(this->enabledWriteModeFiles);
1367      // We only use server-level 1 (for now) to read data
1368      if (serviceType_==CServicesManager::CLIENT)  sendEnabledFiles(this->enabledReadModeFiles);
1369
1370      // Then, send all enabled fields     
1371      sendEnabledFieldsInFiles(this->enabledWriteModeFiles);
1372     
1373      if (serviceType_==CServicesManager::CLIENT) sendEnabledFieldsInFiles(this->enabledReadModeFiles);
1374
1375      // Then, check whether we have domain_ref, axis_ref or scalar_ref attached to the enabled fields
1376      // If any, so send them to server
1377      sendRefDomainsAxisScalars(this->enabledWriteModeFiles); 
1378     
1379      if (serviceType_==CServicesManager::CLIENT) sendRefDomainsAxisScalars(this->enabledReadModeFiles);       
1380
1381      // Check whether enabled fields have grid_ref, if any, send this info to server
1382      sendRefGrid(this->enabledFiles);
1383      // This code may be useful in the future when we want to seperate completely read and write
1384      // sendRefGrid(this->enabledWriteModeFiles);
1385      // if (!hasServer)
1386      //   sendRefGrid(this->enabledReadModeFiles);
1387     
1388      // A grid of enabled fields composed of several components which must be checked then their
1389      // checked attributes should be sent to server
1390      sendGridComponentEnabledFieldsInFiles(this->enabledFiles); // This code can be seperated in two (one for reading, another for writing)
1391
1392      // We have a xml tree on the server side and now, it should be also processed
1393      sendPostProcessing();
1394       
1395      // Finally, we send information of grid itself to server
1396      sendGridEnabledFieldsInFiles(this->enabledWriteModeFiles);       
1397     
1398      if (serviceType_==CServicesManager::CLIENT) sendGridEnabledFieldsInFiles(this->enabledReadModeFiles);       
1399    }
1400    allProcessed = true;
1401
1402
1403    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1404
1405    // There are some processings that should be done after all of above. For example: check mask or index
1406    this->buildFilterGraphOfEnabledFields();
1407   
1408     if (serviceType_==CServicesManager::CLIENT)
1409    {
1410      buildFilterGraphOfFieldsWithReadAccess();
1411      postProcessFilterGraph(); // For coupling in, modify this later
1412    }
1413   
1414    checkGridEnabledFields();   
1415
1416    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1417    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1418
1419    // Nettoyage de l'arborescence
1420    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1421
1422    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1423    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1424   
1425    CTimer::get("Context : close definition").suspend() ;
1426  }
1427  CATCH_DUMP_ATTR
1428
1429 /*!
1430  * Send context attribute and calendar to file server, it must be done once by context file server
1431  * \param[in] client : context client to send   
1432  */ 
1433  void CContext::sendContextToFileServer(CContextClient* client)
1434  {
1435    if (sendToFileServer_done_.count(client)!=0) return ;
1436    else sendToFileServer_done_.insert(client) ;
1437   
1438    this->sendAllAttributesToServer(client); // Send all attributes of current context to server
1439    CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current cale
1440  }
1441
1442  // ym obsolete now to be removed
1443   void CContext::closeDefinition_old(void)
1444   TRY
1445   {
1446     CTimer::get("Context : close definition").resume() ;
1447   
1448    //
1449    postProcessingGlobalAttributes();
1450
1451    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1452
1453    // There are some processings that should be done after all of above. For example: check mask or index
1454    this->buildFilterGraphOfEnabledFields();
1455   
1456     if (serviceType_==CServicesManager::CLIENT)
1457    {
1458      buildFilterGraphOfFieldsWithReadAccess();
1459      postProcessFilterGraph(); // For coupling in, modify this later
1460    }
1461   
1462    checkGridEnabledFields();   
1463
1464    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1465    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1466
1467    // Nettoyage de l'arborescence
1468    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1469
1470    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1471    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1472   
1473    CTimer::get("Context : close definition").suspend() ;
1474   }
1475   CATCH_DUMP_ATTR
1476
1477   vector<CField*> CContext::findAllEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1478   TRY
1479   {
1480     vector<CField*> fields ;
1481     for (unsigned int i = 0; i < activeFiles.size(); i++)
1482     {
1483        const vector<CField*>&& field=activeFiles[i]->getEnabledFields() ;
1484        fields.insert(fields.end(),field.begin(),field.end());
1485     }
1486     return fields ;
1487   }
1488   CATCH_DUMP_ATTR
1489
1490   vector<CField*> CContext::findAllEnabledFieldsInFileOut(const std::vector<CFile*>& activeFiles)
1491   TRY
1492   {
1493     vector<CField*> fields ;
1494     for(auto file : activeFiles)
1495     {
1496        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1497        for(auto field : fieldList) field->setFileOut(file) ;
1498        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1499     }
1500     return fields ;
1501   }
1502   CATCH_DUMP_ATTR
1503
1504   vector<CField*> CContext::findAllEnabledFieldsInFileIn(const std::vector<CFile*>& activeFiles)
1505   TRY
1506   {
1507     vector<CField*> fields ;
1508     for(auto file : activeFiles)
1509     {
1510        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1511        for(auto field : fieldList) field->setFileIn(file) ;
1512        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1513     }
1514     return fields ;
1515   }
1516   CATCH_DUMP_ATTR
1517
1518   vector<CField*> CContext::findAllEnabledFieldsCouplerOut(const std::vector<CCouplerOut*>& activeCouplerOut)
1519   TRY
1520   {
1521     vector<CField*> fields ;
1522     for (auto couplerOut :activeCouplerOut)
1523     {
1524        const vector<CField*>&& fieldList=couplerOut->getEnabledFields() ;
1525        for(auto field : fieldList) field->setCouplerOut(couplerOut) ;
1526        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1527     }
1528     return fields ;
1529   }
1530   CATCH_DUMP_ATTR
1531
1532   vector<CField*> CContext::findAllEnabledFieldsCouplerIn(const std::vector<CCouplerIn*>& activeCouplerIn)
1533   TRY
1534   {
1535     vector<CField*> fields ;
1536     for (auto couplerIn :activeCouplerIn)
1537     {
1538        const vector<CField*>&& fieldList=couplerIn->getEnabledFields() ;
1539        for(auto field : fieldList) field->setCouplerIn(couplerIn) ;
1540        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1541     }
1542     return fields ;
1543   }
1544   CATCH_DUMP_ATTR
1545
1546
1547
1548   void CContext::readAttributesOfEnabledFieldsInReadModeFiles()
1549   TRY
1550   {
1551      for (unsigned int i = 0; i < this->enabledReadModeFiles.size(); ++i)
1552        (void)this->enabledReadModeFiles[i]->readAttributesOfEnabledFieldsInReadMode();
1553   }
1554   CATCH_DUMP_ATTR
1555
1556   void CContext::sendGridComponentEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1557   TRY
1558   {
1559     int size = activeFiles.size();
1560     for (int i = 0; i < size; ++i)
1561     {       
1562       activeFiles[i]->sendGridComponentOfEnabledFields();
1563     }
1564   }
1565   CATCH_DUMP_ATTR
1566
1567   /*!
1568      Send active (enabled) fields in file from a client to others
1569      \param [in] activeFiles files contains enabled fields to send
1570   */
1571   void CContext::sendGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1572   TRY
1573   {
1574     int size = activeFiles.size();
1575     for (int i = 0; i < size; ++i)
1576     {       
1577       activeFiles[i]->sendGridOfEnabledFields();
1578     }
1579   }
1580   CATCH_DUMP_ATTR
1581
1582   void CContext::checkGridEnabledFields()
1583   TRY
1584   {
1585     int size = enabledFiles.size();
1586     for (int i = 0; i < size; ++i)
1587     {
1588       enabledFiles[i]->checkGridOfEnabledFields();       
1589     }
1590
1591     size = enabledCouplerOut.size();
1592     for (int i = 0; i < size; ++i)
1593     {
1594       enabledCouplerOut[i]->checkGridOfEnabledFields();       
1595     }
1596   }
1597   CATCH_DUMP_ATTR
1598
1599   /*!
1600      Check grid of active (enabled) fields in file
1601      \param [in] activeFiles files contains enabled fields whose grid needs checking
1602   */
1603   void CContext::checkGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1604   TRY
1605   {
1606     int size = activeFiles.size();
1607     for (int i = 0; i < size; ++i)
1608     {
1609       activeFiles[i]->checkGridOfEnabledFields();       
1610     }
1611   }
1612   CATCH_DUMP_ATTR
1613
1614    /*!
1615      Go up the hierachical tree via field_ref and do check of attributes of fields
1616      This can be done in a client then all computed information will be sent from this client to others
1617      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1618   */
1619   void CContext::solveOnlyRefOfEnabledFields(void)
1620   TRY
1621   {
1622     int size = this->enabledFiles.size();
1623     for (int i = 0; i < size; ++i)
1624     {
1625       this->enabledFiles[i]->solveOnlyRefOfEnabledFields();
1626     }
1627
1628     for (int i = 0; i < size; ++i)
1629     {
1630       this->enabledFiles[i]->generateNewTransformationGridDest();
1631     }
1632
1633     size = this->enabledCouplerOut.size();
1634     for (int i = 0; i < size; ++i)
1635     {
1636       this->enabledCouplerOut[i]->solveOnlyRefOfEnabledFields();
1637     }
1638
1639     for (int i = 0; i < size; ++i)
1640     {
1641       this->enabledCouplerOut[i]->generateNewTransformationGridDest();
1642     }
1643   }
1644   CATCH_DUMP_ATTR
1645
1646    /*!
1647      Go up the hierachical tree via field_ref and do check of attributes of fields.
1648      The transformation can be done in this step.
1649      All computed information will be sent from this client to others.
1650      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1651   */
1652   void CContext::solveAllRefOfEnabledFieldsAndTransform(void)
1653   TRY
1654   {
1655     int size = this->enabledFiles.size();
1656     for (int i = 0; i < size; ++i)
1657     {
1658       this->enabledFiles[i]->solveAllRefOfEnabledFieldsAndTransform();
1659     }
1660
1661     size = this->enabledCouplerOut.size();
1662     for (int i = 0; i < size; ++i)
1663     {
1664       this->enabledCouplerOut[i]->solveAllRefOfEnabledFieldsAndTransform();
1665     }
1666
1667   }
1668   CATCH_DUMP_ATTR
1669
1670   void CContext::buildFilterGraphOfEnabledFields()
1671   TRY
1672   {
1673     int size = this->enabledFiles.size();
1674     for (int i = 0; i < size; ++i)
1675     {
1676       this->enabledFiles[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1677     }
1678
1679     size = this->enabledCouplerOut.size();
1680     for (int i = 0; i < size; ++i)
1681     {
1682       this->enabledCouplerOut[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1683     }
1684   }
1685   CATCH_DUMP_ATTR
1686
1687   void CContext::postProcessFilterGraph()
1688   TRY
1689   {
1690     int size = enabledFiles.size();
1691     for (int i = 0; i < size; ++i)
1692     {
1693        enabledFiles[i]->postProcessFilterGraph();
1694     }
1695   }
1696   CATCH_DUMP_ATTR
1697
1698   void CContext::startPrefetchingOfEnabledReadModeFiles()
1699   TRY
1700   {
1701     int size = enabledReadModeFiles.size();
1702     for (int i = 0; i < size; ++i)
1703     {
1704        enabledReadModeFiles[i]->prefetchEnabledReadModeFields();
1705     }
1706   }
1707   CATCH_DUMP_ATTR
1708
1709   void CContext::doPreTimestepOperationsForEnabledReadModeFiles()
1710   TRY
1711   {
1712     int size = enabledReadModeFiles.size();
1713     for (int i = 0; i < size; ++i)
1714     {
1715        enabledReadModeFiles[i]->doPreTimestepOperationsForEnabledReadModeFields();
1716     }
1717   }
1718   CATCH_DUMP_ATTR
1719
1720   void CContext::doPostTimestepOperationsForEnabledReadModeFiles()
1721   TRY
1722   {
1723     int size = enabledReadModeFiles.size();
1724     for (int i = 0; i < size; ++i)
1725     {
1726        enabledReadModeFiles[i]->doPostTimestepOperationsForEnabledReadModeFields();
1727     }
1728   }
1729   CATCH_DUMP_ATTR
1730
1731  void CContext::findFieldsWithReadAccess(void)
1732  TRY
1733  {
1734    fieldsWithReadAccess_.clear();
1735    const vector<CField*> allFields = CField::getAll();
1736    for (size_t i = 0; i < allFields.size(); ++i)
1737    {
1738      CField* field = allFields[i];
1739      if (!field->read_access.isEmpty() && field->read_access && (field->enabled.isEmpty() || field->enabled))
1740      {
1741        fieldsWithReadAccess_.push_back(field);
1742        field->setModelOut() ;
1743      }
1744    }
1745  }
1746  CATCH_DUMP_ATTR
1747
1748  void CContext::solveAllRefOfFieldsWithReadAccess()
1749  TRY
1750  {
1751    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1752      fieldsWithReadAccess_[i]->solveAllReferenceEnabledField(false);
1753  }
1754  CATCH_DUMP_ATTR
1755
1756  void CContext::buildFilterGraphOfFieldsWithReadAccess()
1757  TRY
1758  {
1759    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1760      fieldsWithReadAccess_[i]->buildFilterGraph(garbageCollector, true);
1761  }
1762  CATCH_DUMP_ATTR
1763
1764   void CContext::solveAllInheritance(bool apply)
1765   TRY
1766   {
1767     // Résolution des héritages descendants (càd des héritages de groupes)
1768     // pour chacun des contextes.
1769      solveDescInheritance(apply);
1770
1771     // Résolution des héritages par référence au niveau des fichiers.
1772      const vector<CFile*> allFiles=CFile::getAll();
1773      const vector<CCouplerIn*> allCouplerIn=CCouplerIn::getAll();
1774      const vector<CCouplerOut*> allCouplerOut=CCouplerOut::getAll();
1775      const vector<CGrid*> allGrids= CGrid::getAll();
1776
1777      if (serviceType_==CServicesManager::CLIENT)
1778      {
1779        for (unsigned int i = 0; i < allFiles.size(); i++)
1780          allFiles[i]->solveFieldRefInheritance(apply);
1781
1782        for (unsigned int i = 0; i < allCouplerIn.size(); i++)
1783          allCouplerIn[i]->solveFieldRefInheritance(apply);
1784
1785        for (unsigned int i = 0; i < allCouplerOut.size(); i++)
1786          allCouplerOut[i]->solveFieldRefInheritance(apply);
1787      }
1788
1789      unsigned int vecSize = allGrids.size();
1790      unsigned int i = 0;
1791      for (i = 0; i < vecSize; ++i)
1792        allGrids[i]->solveElementsRefInheritance(apply);
1793
1794   }
1795  CATCH_DUMP_ATTR
1796
1797   void CContext::findEnabledFiles(void)
1798   TRY
1799   {
1800      const std::vector<CFile*> allFiles = CFile::getAll();
1801      const CDate& initDate = calendar->getInitDate();
1802
1803      for (unsigned int i = 0; i < allFiles.size(); i++)
1804         if (!allFiles[i]->enabled.isEmpty()) // Si l'attribut 'enabled' est défini.
1805         {
1806            if (allFiles[i]->enabled.getValue()) // Si l'attribut 'enabled' est fixé à vrai.
1807            {
1808              if (allFiles[i]->output_freq.isEmpty())
1809              {
1810                 ERROR("CContext::findEnabledFiles()",
1811                     << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1812                     <<" \".")
1813              }
1814              if ((initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1815              {
1816                error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1817                    << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1818                    <<"\" is less than the time step. File will not be written."<<endl;
1819              }
1820              else
1821               enabledFiles.push_back(allFiles[i]);
1822            }
1823         }
1824         else
1825         {
1826           if (allFiles[i]->output_freq.isEmpty())
1827           {
1828              ERROR("CContext::findEnabledFiles()",
1829                  << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1830                  <<" \".")
1831           }
1832           if ( (initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1833           {
1834             error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1835                 << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1836                 <<"\" is less than the time step. File will not be written."<<endl;
1837           }
1838           else
1839             enabledFiles.push_back(allFiles[i]); // otherwise true by default
1840         }
1841
1842      if (enabledFiles.size() == 0)
1843         DEBUG(<<"Aucun fichier ne va être sorti dans le contexte nommé \""
1844               << getId() << "\" !");
1845
1846   }
1847   CATCH_DUMP_ATTR
1848
1849   void CContext::findEnabledCouplerIn(void)
1850   TRY
1851   {
1852      const std::vector<CCouplerIn*> allCouplerIn = CCouplerIn::getAll();
1853      bool enabled ;
1854      for (size_t i = 0; i < allCouplerIn.size(); i++)
1855      {
1856        if (allCouplerIn[i]->enabled.isEmpty()) enabled=true ;
1857        else enabled=allCouplerIn[i]->enabled ;
1858        if (enabled) enabledCouplerIn.push_back(allCouplerIn[i]) ;
1859      }
1860   }
1861   CATCH_DUMP_ATTR
1862
1863   void CContext::findEnabledCouplerOut(void)
1864   TRY
1865   {
1866      const std::vector<CCouplerOut*> allCouplerOut = CCouplerOut::getAll();
1867      bool enabled ;
1868      for (size_t i = 0; i < allCouplerOut.size(); i++)
1869      {
1870        if (allCouplerOut[i]->enabled.isEmpty()) enabled=true ;
1871        else enabled=allCouplerOut[i]->enabled ;
1872        if (enabled) enabledCouplerOut.push_back(allCouplerOut[i]) ;
1873      }
1874   }
1875   CATCH_DUMP_ATTR
1876
1877
1878
1879
1880   void CContext::distributeFiles(const vector<CFile*>& files)
1881   TRY
1882   {
1883     bool distFileMemory=false ;
1884     distFileMemory=CXios::getin<bool>("server2_dist_file_memory", distFileMemory);
1885
1886     if (distFileMemory) distributeFileOverMemoryBandwith(files) ;
1887     else distributeFileOverBandwith(files) ;
1888   }
1889   CATCH_DUMP_ATTR
1890
1891   void CContext::distributeFileOverBandwith(const vector<CFile*>& files)
1892   TRY
1893   {
1894     double eps=std::numeric_limits<double>::epsilon()*10 ;
1895     
1896     std::ofstream ofs(("distribute_file_"+getId()+".dat").c_str(), std::ofstream::out);
1897     int nbPools = clientPrimServer.size();
1898
1899     // (1) Find all enabled files in write mode
1900     // for (int i = 0; i < this->enabledFiles.size(); ++i)
1901     // {
1902     //   if (enabledFiles[i]->mode.isEmpty() || (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
1903     //    enabledWriteModeFiles.push_back(enabledFiles[i]);
1904     // }
1905
1906     // (2) Estimate the data volume for each file
1907     int size = files.size();
1908     std::vector<std::pair<double, CFile*> > dataSizeMap;
1909     double dataPerPool = 0;
1910     int nfield=0 ;
1911     ofs<<size<<endl ;
1912     for (size_t i = 0; i < size; ++i)
1913     {
1914       CFile* file = files[i];
1915       ofs<<file->getId()<<endl ;
1916       StdSize dataSize=0;
1917       std::vector<CField*> enabledFields = file->getEnabledFields();
1918       size_t numEnabledFields = enabledFields.size();
1919       ofs<<numEnabledFields<<endl ;
1920       for (size_t j = 0; j < numEnabledFields; ++j)
1921       {
1922         dataSize += enabledFields[j]->getGlobalWrittenSize() ;
1923         ofs<<enabledFields[j]->getGrid()->getId()<<endl ;
1924         ofs<<enabledFields[j]->getGlobalWrittenSize()<<endl ;
1925       }
1926       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
1927       double dataSizeSec= dataSize/ outFreqSec;
1928       ofs<<dataSizeSec<<endl ;
1929       nfield++ ;
1930// add epsilon*nField to dataSizeSec in order to  preserve reproductive ordering when sorting
1931       dataSizeMap.push_back(make_pair(dataSizeSec + dataSizeSec * eps * nfield , file));
1932       dataPerPool += dataSizeSec;
1933     }
1934     dataPerPool /= nbPools;
1935     std::sort(dataSizeMap.begin(), dataSizeMap.end());
1936
1937     // (3) Assign contextClient to each enabled file
1938
1939     std::multimap<double,int> poolDataSize ;
1940// multimap is not garanty to preserve stable sorting in c++98 but it seems it does for c++11
1941
1942     int j;
1943     double dataSize ;
1944     for (j = 0 ; j < nbPools ; ++j) poolDataSize.insert(std::pair<double,int>(0.,j)) ; 
1945             
1946     for (int i = dataSizeMap.size()-1; i >= 0; --i)
1947     {
1948       dataSize=(*poolDataSize.begin()).first ;
1949       j=(*poolDataSize.begin()).second ;
1950       dataSizeMap[i].second->setContextClient(clientPrimServer[j]);
1951       dataSize+=dataSizeMap[i].first;
1952       poolDataSize.erase(poolDataSize.begin()) ;
1953       poolDataSize.insert(std::pair<double,int>(dataSize,j)) ; 
1954     }
1955
1956     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 ;
1957   }
1958   CATCH_DUMP_ATTR
1959
1960   void CContext::distributeFileOverMemoryBandwith(const vector<CFile*>& filesList)
1961   TRY
1962   {
1963     int nbPools = clientPrimServer.size();
1964     double ratio=0.5 ;
1965     ratio=CXios::getin<double>("server2_dist_file_memory_ratio", ratio);
1966
1967     int nFiles = filesList.size();
1968     vector<SDistFile> files(nFiles);
1969     vector<SDistGrid> grids;
1970     map<string,int> gridMap ;
1971     string gridId; 
1972     int gridIndex=0 ;
1973
1974     for (size_t i = 0; i < nFiles; ++i)
1975     {
1976       StdSize dataSize=0;
1977       CFile* file = filesList[i];
1978       std::vector<CField*> enabledFields = file->getEnabledFields();
1979       size_t numEnabledFields = enabledFields.size();
1980
1981       files[i].id_=file->getId() ;
1982       files[i].nbGrids_=numEnabledFields;
1983       files[i].assignedGrid_ = new int[files[i].nbGrids_] ;
1984         
1985       for (size_t j = 0; j < numEnabledFields; ++j)
1986       {
1987         gridId=enabledFields[j]->getGrid()->getId() ;
1988         if (gridMap.find(gridId)==gridMap.end())
1989         {
1990            gridMap[gridId]=gridIndex  ;
1991            SDistGrid newGrid; 
1992            grids.push_back(newGrid) ;
1993            gridIndex++ ;
1994         }
1995         files[i].assignedGrid_[j]=gridMap[gridId] ;
1996         grids[files[i].assignedGrid_[j]].size_=enabledFields[j]->getGlobalWrittenSize() ;
1997         dataSize += enabledFields[j]->getGlobalWrittenSize() ; // usefull
1998       }
1999       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
2000       files[i].bandwith_= dataSize/ outFreqSec ;
2001     }
2002
2003     double bandwith=0 ;
2004     double memory=0 ;
2005   
2006     for(int i=0; i<nFiles; i++)  bandwith+=files[i].bandwith_ ;
2007     for(int i=0; i<nFiles; i++)  files[i].bandwith_ = files[i].bandwith_/bandwith * ratio ;
2008
2009     for(int i=0; i<grids.size(); i++)  memory+=grids[i].size_ ;
2010     for(int i=0; i<grids.size(); i++)  grids[i].size_ = grids[i].size_ / memory * (1.0-ratio) ;
2011       
2012     distributeFileOverServer2(nbPools, grids.size(), &grids[0], nFiles, &files[0]) ;
2013
2014     vector<double> memorySize(nbPools,0.) ;
2015     vector< set<int> > serverGrids(nbPools) ;
2016     vector<double> bandwithSize(nbPools,0.) ;
2017       
2018     for (size_t i = 0; i < nFiles; ++i)
2019     {
2020       bandwithSize[files[i].assignedServer_] += files[i].bandwith_* bandwith /ratio ;
2021       for(int j=0 ; j<files[i].nbGrids_;j++)
2022       {
2023         if (serverGrids[files[i].assignedServer_].find(files[i].assignedGrid_[j]) == serverGrids[files[i].assignedServer_].end())
2024         {
2025           memorySize[files[i].assignedServer_]+= grids[files[i].assignedGrid_[j]].size_ * memory / (1.0-ratio);
2026           serverGrids[files[i].assignedServer_].insert(files[i].assignedGrid_[j]) ;
2027         }
2028       }
2029       filesList[i]->setContextClient(clientPrimServer[files[i].assignedServer_]) ;
2030       delete [] files[i].assignedGrid_ ;
2031     }
2032
2033     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 ;
2034     for (int i = 0; i < nbPools; ++i) info(100)<<"Pool server level2 "<<i<<"   assigned grid memory "<<memorySize[i]*100/1024./1024.<<" Mb"<<endl ;
2035
2036   }
2037   CATCH_DUMP_ATTR
2038
2039   /*!
2040      Find all files in write mode
2041   */
2042   void CContext::findEnabledWriteModeFiles(void)
2043   TRY
2044   {
2045     int size = this->enabledFiles.size();
2046     for (int i = 0; i < size; ++i)
2047     {
2048       if (enabledFiles[i]->mode.isEmpty() || 
2049          (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
2050        enabledWriteModeFiles.push_back(enabledFiles[i]);
2051     }
2052   }
2053   CATCH_DUMP_ATTR
2054
2055   /*!
2056      Find all files in read mode
2057   */
2058   void CContext::findEnabledReadModeFiles(void)
2059   TRY
2060   {
2061     int size = this->enabledFiles.size();
2062     for (int i = 0; i < size; ++i)
2063     {
2064       if (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::read)
2065        enabledReadModeFiles.push_back(enabledFiles[i]);
2066     }
2067   }
2068   CATCH_DUMP_ATTR
2069
2070   void CContext::closeAllFile(void)
2071   TRY
2072   {
2073     std::vector<CFile*>::const_iterator
2074            it = this->enabledFiles.begin(), end = this->enabledFiles.end();
2075
2076     for (; it != end; it++)
2077     {
2078       info(30)<<"Closing File : "<<(*it)->getId()<<endl;
2079       (*it)->close();
2080     }
2081   }
2082   CATCH_DUMP_ATTR
2083
2084   /*!
2085   \brief Dispatch event received from client
2086      Whenever a message is received in buffer of server, it will be processed depending on
2087   its event type. A new event type should be added in the switch list to make sure
2088   it processed on server side.
2089   \param [in] event: Received message
2090   */
2091   bool CContext::dispatchEvent(CEventServer& event)
2092   TRY
2093   {
2094
2095      if (SuperClass::dispatchEvent(event)) return true;
2096      else
2097      {
2098        switch(event.type)
2099        {
2100           case EVENT_ID_CLOSE_DEFINITION :
2101             recvCloseDefinition(event);
2102             return true;
2103             break;
2104           case EVENT_ID_UPDATE_CALENDAR:
2105             recvUpdateCalendar(event);
2106             return true;
2107             break;
2108           case EVENT_ID_CREATE_FILE_HEADER :
2109             recvCreateFileHeader(event);
2110             return true;
2111             break;
2112           case EVENT_ID_POST_PROCESS:
2113             recvPostProcessing(event);
2114             return true;
2115           case EVENT_ID_SEND_REGISTRY:
2116             recvRegistry(event);
2117             return true;
2118             break;
2119           case EVENT_ID_POST_PROCESS_GLOBAL_ATTRIBUTES:
2120             recvPostProcessingGlobalAttributes(event);
2121             return true;
2122             break;
2123           case EVENT_ID_PROCESS_GRID_ENABLED_FIELDS:
2124             recvProcessingGridOfEnabledFields(event);
2125             return true;
2126             break;
2127           case EVENT_ID_COUPLER_IN_READY:
2128             recvCouplerInReady(event);
2129             return true;
2130             break;
2131           case EVENT_ID_COUPLER_IN_CLOSE_DEFINITION:
2132             recvCouplerInCloseDefinition(event);
2133             return true;
2134             break;
2135           case EVENT_ID_COUPLER_IN_CONTEXT_FINALIZED:
2136             recvCouplerInContextFinalized(event);
2137             return true;
2138             break; 
2139           default :
2140             ERROR("bool CContext::dispatchEvent(CEventServer& event)",
2141                    <<"Unknown Event");
2142           return false;
2143         }
2144      }
2145   }
2146   CATCH
2147
2148   //! Client side: Send a message to server to make it close
2149   // ym obsolete
2150   void CContext::sendCloseDefinition(void)
2151   TRY
2152   {
2153    int nbSrvPools ;
2154    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2155    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2156    else nbSrvPools = 0 ;
2157    CContextClient* contextClientTmp ;
2158
2159    for (int i = 0; i < nbSrvPools; ++i)
2160     {
2161       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2162       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2163       CEventClient event(getType(),EVENT_ID_CLOSE_DEFINITION);
2164       if (contextClientTmp->isServerLeader())
2165       {
2166         CMessage msg;
2167         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2168         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2169           event.push(*itRank,1,msg);
2170         contextClientTmp->sendEvent(event);
2171       }
2172       else contextClientTmp->sendEvent(event);
2173     }
2174   }
2175   CATCH_DUMP_ATTR
2176   
2177   //  ! Client side: Send a message to server to make it close
2178   void CContext::sendCloseDefinition(CContextClient* client)
2179   TRY
2180   {
2181      if (sendCloseDefinition_done_.count(client)!=0) return ;
2182      else sendCloseDefinition_done_.insert(client) ;
2183
2184      CEventClient event(getType(),EVENT_ID_CLOSE_DEFINITION);
2185      if (client->isServerLeader())
2186      {
2187        CMessage msg;
2188        for(auto rank : client->getRanksServerLeader()) event.push(rank,1,msg);
2189        client->sendEvent(event);
2190      }
2191     else client->sendEvent(event);
2192   }
2193   CATCH_DUMP_ATTR
2194
2195   //! Server side: Receive a message of client announcing a context close
2196   void CContext::recvCloseDefinition(CEventServer& event)
2197   TRY
2198   {
2199      CBufferIn* buffer=event.subEvents.begin()->buffer;
2200      getCurrent()->closeDefinition();
2201   }
2202   CATCH
2203
2204   //! Client side: Send a message to update calendar in each time step
2205   void CContext::sendUpdateCalendar(int step)
2206   TRY
2207   {
2208     CEventClient event(getType(),EVENT_ID_UPDATE_CALENDAR);
2209     for(auto client : slaveServers_) 
2210     {
2211       if (client->isServerLeader())
2212       {
2213         CMessage msg;
2214         msg<<step;
2215         for (auto& rank : client->getRanksServerLeader() ) event.push(rank,1,msg);
2216         client->sendEvent(event);
2217       }
2218       else client->sendEvent(event);
2219     }
2220   }
2221   CATCH_DUMP_ATTR
2222
2223   //! Server side: Receive a message of client annoucing calendar update
2224   void CContext::recvUpdateCalendar(CEventServer& event)
2225   TRY
2226   {
2227      CBufferIn* buffer=event.subEvents.begin()->buffer;
2228      getCurrent()->recvUpdateCalendar(*buffer);
2229   }
2230   CATCH
2231
2232   //! Server side: Receive a message of client annoucing calendar update
2233   void CContext::recvUpdateCalendar(CBufferIn& buffer)
2234   TRY
2235   {
2236      int step;
2237      buffer>>step;
2238      updateCalendar(step);
2239      if (serviceType_==CServicesManager::GATHERER)
2240      {       
2241        sendUpdateCalendar(step);
2242      }
2243   }
2244   CATCH_DUMP_ATTR
2245
2246   //! Client side: Send a message to create header part of netcdf file
2247   void CContext::sendCreateFileHeader(void)
2248   TRY
2249   {
2250     int nbSrvPools ;
2251     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2252     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2253     else nbSrvPools = 0 ;
2254     CContextClient* contextClientTmp ;
2255
2256     for (int i = 0; i < nbSrvPools; ++i)
2257     {
2258       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2259       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2260       CEventClient event(getType(),EVENT_ID_CREATE_FILE_HEADER);
2261
2262       if (contextClientTmp->isServerLeader())
2263       {
2264         CMessage msg;
2265         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2266         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2267           event.push(*itRank,1,msg) ;
2268         contextClientTmp->sendEvent(event);
2269       }
2270       else contextClientTmp->sendEvent(event);
2271     }
2272   }
2273   CATCH_DUMP_ATTR
2274
2275   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2276   void CContext::recvCreateFileHeader(CEventServer& event)
2277   TRY
2278   {
2279      CBufferIn* buffer=event.subEvents.begin()->buffer;
2280      getCurrent()->recvCreateFileHeader(*buffer);
2281   }
2282   CATCH
2283
2284   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2285   void CContext::recvCreateFileHeader(CBufferIn& buffer)
2286   TRY
2287   {
2288      if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER) 
2289        createFileHeader();
2290   }
2291   CATCH_DUMP_ATTR
2292
2293   //! Client side: Send a message to do some post processing on server
2294   void CContext::sendProcessingGridOfEnabledFields()
2295   TRY
2296   {
2297     int nbSrvPools ;
2298     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2299     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2300     else nbSrvPools = 0 ;
2301     CContextClient* contextClientTmp ;
2302
2303     for (int i = 0; i < nbSrvPools; ++i)
2304     {
2305       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2306       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2307
2308       CEventClient event(getType(),EVENT_ID_PROCESS_GRID_ENABLED_FIELDS);
2309
2310       if (contextClientTmp->isServerLeader())
2311       {
2312         CMessage msg;
2313         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2314         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2315           event.push(*itRank,1,msg);
2316         contextClientTmp->sendEvent(event);
2317       }
2318       else contextClientTmp->sendEvent(event);
2319     }
2320   }
2321   CATCH_DUMP_ATTR
2322
2323   //! Server side: Receive a message to do some post processing
2324   void CContext::recvProcessingGridOfEnabledFields(CEventServer& event)
2325   TRY
2326   {
2327      CBufferIn* buffer=event.subEvents.begin()->buffer;
2328      // nothing to do, no call ??!!   
2329   }
2330   CATCH
2331
2332   //! Client side: Send a message to do some post processing on server
2333   void CContext::sendPostProcessing()
2334   TRY
2335   {
2336     int nbSrvPools ;
2337     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2338     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2339     else nbSrvPools = 0 ;
2340     CContextClient* contextClientTmp ;
2341
2342     for (int i = 0; i < nbSrvPools; ++i)
2343     {
2344       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2345       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2346       CEventClient event(getType(),EVENT_ID_POST_PROCESS);
2347       if (contextClientTmp->isServerLeader())
2348       {
2349         CMessage msg;
2350         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2351         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2352         event.push(*itRank,1,msg);
2353         contextClientTmp->sendEvent(event);
2354       }
2355       else contextClientTmp->sendEvent(event);
2356     }
2357   }
2358   CATCH_DUMP_ATTR
2359
2360   //! Server side: Receive a message to do some post processing
2361   void CContext::recvPostProcessing(CEventServer& event)
2362   TRY
2363   {
2364      CBufferIn* buffer=event.subEvents.begin()->buffer;
2365      getCurrent()->recvPostProcessing(*buffer);
2366   }
2367   CATCH
2368
2369   //! Server side: Receive a message to do some post processing
2370   void CContext::recvPostProcessing(CBufferIn& buffer)
2371   TRY
2372   {
2373      CCalendarWrapper::get(CCalendarWrapper::GetDefName())->createCalendar();
2374      postProcessing();
2375   }
2376   CATCH_DUMP_ATTR
2377
2378 
2379   /*!
2380   \brief Do some simple post processings after parsing xml file
2381      After the xml file (iodef.xml) is parsed, it is necessary to build all relations among
2382   created object, e.g: inhertance among fields, domain, axis. After that, all fiels as well as their parents (reference fields),
2383   which will be written out into netcdf files, are processed
2384   */
2385   void CContext::postProcessing()
2386   TRY
2387   {
2388     if (isPostProcessed) return;
2389
2390      // Make sure the calendar was correctly created
2391      if (!calendar)
2392        ERROR("CContext::postProcessing()", << "A calendar must be defined for the context \"" << getId() << "!\"")
2393      else if (calendar->getTimeStep() == NoneDu)
2394        ERROR("CContext::postProcessing()", << "A timestep must be defined for the context \"" << getId() << "!\"")
2395      // Calendar first update to set the current date equals to the start date
2396      calendar->update(0);
2397
2398      // Find all inheritance in xml structure
2399      this->solveAllInheritance();
2400
2401//      ShowTree(info(10));
2402
2403      // Check if some axis, domains or grids are eligible to for compressed indexed output.
2404      // Warning: This must be done after solving the inheritance and before the rest of post-processing
2405      checkAxisDomainsGridsEligibilityForCompressedOutput();      // only for field written on IO_SERVER service ????
2406
2407      // Check if some automatic time series should be generated
2408      // Warning: This must be done after solving the inheritance and before the rest of post-processing     
2409      prepareTimeseries();
2410
2411      //Initialisation du vecteur 'enabledFiles' contenant la liste des fichiers à sortir.
2412      findEnabledFiles();
2413      findEnabledWriteModeFiles();
2414      findEnabledReadModeFiles();
2415      findEnabledCouplerIn();
2416      findEnabledCouplerOut();
2417      createCouplerInterCommunicator() ;
2418
2419      // Find all enabled fields of each file     
2420      const vector<CField*>&& fileOutField = findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
2421      const vector<CField*>&& fileInField = findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
2422      const vector<CField*>&& CouplerOutField = findAllEnabledFieldsCouplerOut(this->enabledCouplerOut);
2423      const vector<CField*>&& CouplerInField = findAllEnabledFieldsCouplerIn(this->enabledCouplerIn);
2424
2425
2426
2427      // For now, only read files with client and only one level server
2428      // if (hasClient && !hasServer) findEnabledReadModeFiles();     
2429
2430
2431      // For now, only read files with client and only one level server
2432      // if (hasClient && !hasServer)
2433      //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
2434
2435      if (serviceType_==CServicesManager::CLIENT)
2436      {
2437        initReadFiles();
2438        // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
2439        this->readAttributesOfEnabledFieldsInReadModeFiles();
2440      }
2441
2442      // Only search and rebuild all reference objects of enable fields, don't transform
2443      this->solveOnlyRefOfEnabledFields();
2444
2445      // Search and rebuild all reference object of enabled fields, and transform
2446      this->solveAllRefOfEnabledFieldsAndTransform();
2447
2448      // Find all fields with read access from the public API
2449      if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
2450      // and solve the all reference for them
2451      if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
2452
2453      isPostProcessed = true;
2454   }
2455   CATCH_DUMP_ATTR
2456
2457   void CContext::createCouplerInterCommunicator(void)
2458   TRY
2459   {
2460      int rank=this->getIntraCommRank() ;
2461      map<string,list<CCouplerOut*>> listCouplerOut ; 
2462      map<string,list<CCouplerIn*>> listCouplerIn ; 
2463
2464      for(auto couplerOut : enabledCouplerOut) listCouplerOut[couplerOut->getCouplingContextId()].push_back(couplerOut) ;
2465      for(auto couplerIn : enabledCouplerIn) listCouplerIn[couplerIn->getCouplingContextId()].push_back(couplerIn) ;
2466
2467      CCouplerManager* couplerManager = CXios::getCouplerManager() ;
2468      if (rank==0)
2469      {
2470        for(auto couplerOut : listCouplerOut) couplerManager->registerCoupling(this->getContextId(),couplerOut.first) ;
2471        for(auto couplerIn : listCouplerIn) couplerManager->registerCoupling(couplerIn.first,this->getContextId()) ;
2472      }
2473
2474      do
2475      {
2476        for(auto couplerOut : listCouplerOut) 
2477        {
2478          bool isNextCoupling ;
2479          if (rank==0) isNextCoupling = couplerManager->isNextCoupling(this->getContextId(),couplerOut.first) ;
2480          MPI_Bcast(&isNextCoupling,1,MPI_C_BOOL, 0, getIntraComm()) ; 
2481          if (isNextCoupling) 
2482          {
2483            addCouplingChanel(couplerOut.first, true) ;
2484            listCouplerOut.erase(couplerOut.first) ;
2485            break ;
2486          }           
2487        }
2488        for(auto couplerIn : listCouplerIn) 
2489        {
2490          bool isNextCoupling ;
2491          if (rank==0) isNextCoupling = couplerManager->isNextCoupling(couplerIn.first,this->getContextId());
2492          MPI_Bcast(&isNextCoupling,1,MPI_C_BOOL, 0, getIntraComm()) ; 
2493          if (isNextCoupling) 
2494          {
2495            addCouplingChanel(couplerIn.first, false) ;
2496            listCouplerIn.erase(couplerIn.first) ;
2497            break ;
2498          }           
2499        }
2500
2501      } while (!listCouplerOut.empty() || !listCouplerIn.empty()) ;
2502
2503   }
2504   CATCH_DUMP_ATTR
2505
2506 
2507     //! Client side: Send infomation of active files (files are enabled to write out)
2508   void CContext::sendEnabledFiles(const std::vector<CFile*>& activeFiles)
2509   TRY
2510   {
2511     int size = activeFiles.size();
2512
2513     // In a context, each type has a root definition, e.g: axis, domain, field.
2514     // Every object must be a child of one of these root definition. In this case
2515     // all new file objects created on server must be children of the root "file_definition"
2516     StdString fileDefRoot("file_definition");
2517     CFileGroup* cfgrpPtr = CFileGroup::get(fileDefRoot);
2518
2519     for (int i = 0; i < size; ++i)
2520     {
2521       CFile* f = activeFiles[i];
2522       cfgrpPtr->sendCreateChild(f->getId(),f->getContextClient());
2523       f->sendAllAttributesToServer(f->getContextClient());
2524       f->sendAddAllVariables(f->getContextClient());
2525     }
2526   }
2527   CATCH_DUMP_ATTR
2528
2529   //! Client side: Send information of active fields (ones are written onto files)
2530   void CContext::sendEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
2531   TRY
2532   {
2533     int size = activeFiles.size();
2534     for (int i = 0; i < size; ++i)
2535     {
2536       activeFiles[i]->sendEnabledFields(activeFiles[i]->getContextClient());
2537     }
2538   }
2539   CATCH_DUMP_ATTR
2540
2541   //! Client side: Check if the defined axis, domains and grids are eligible for compressed indexed output
2542   void CContext::checkAxisDomainsGridsEligibilityForCompressedOutput()
2543   TRY
2544   {
2545     if (!(serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)) return;
2546
2547     const vector<CAxis*> allAxis = CAxis::getAll();
2548     for (vector<CAxis*>::const_iterator it = allAxis.begin(); it != allAxis.end(); it++)
2549       (*it)->checkEligibilityForCompressedOutput();
2550
2551     const vector<CDomain*> allDomains = CDomain::getAll();
2552     for (vector<CDomain*>::const_iterator it = allDomains.begin(); it != allDomains.end(); it++)
2553       (*it)->checkEligibilityForCompressedOutput();
2554
2555     const vector<CGrid*> allGrids = CGrid::getAll();
2556     for (vector<CGrid*>::const_iterator it = allGrids.begin(); it != allGrids.end(); it++)
2557       (*it)->checkEligibilityForCompressedOutput();
2558   }
2559   CATCH_DUMP_ATTR
2560
2561   //! Client side: Prepare the timeseries by adding the necessary files
2562   void CContext::prepareTimeseries()
2563   TRY
2564   {
2565     const std::vector<CFile*> allFiles = CFile::getAll();
2566     for (size_t i = 0; i < allFiles.size(); i++)
2567     {
2568       CFile* file = allFiles[i];
2569
2570       std::vector<CVariable*> fileVars, fieldVars, vars = file->getAllVariables();
2571       for (size_t k = 0; k < vars.size(); k++)
2572       {
2573         CVariable* var = vars[k];
2574
2575         if (var->ts_target.isEmpty()
2576              || var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both)
2577           fileVars.push_back(var);
2578
2579         if (!var->ts_target.isEmpty()
2580              && (var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both))
2581           fieldVars.push_back(var);
2582       }
2583
2584       if (!file->timeseries.isEmpty() && file->timeseries != CFile::timeseries_attr::none)
2585       {
2586         StdString fileNameStr("%file_name%") ;
2587         StdString tsPrefix = !file->ts_prefix.isEmpty() ? file->ts_prefix : fileNameStr ;
2588         
2589         StdString fileName=file->getFileOutputName();
2590         size_t pos=tsPrefix.find(fileNameStr) ;
2591         while (pos!=std::string::npos)
2592         {
2593           tsPrefix=tsPrefix.replace(pos,fileNameStr.size(),fileName) ;
2594           pos=tsPrefix.find(fileNameStr) ;
2595         }
2596       
2597         const std::vector<CField*> allFields = file->getAllFields();
2598         for (size_t j = 0; j < allFields.size(); j++)
2599         {
2600           CField* field = allFields[j];
2601
2602           if (!field->ts_enabled.isEmpty() && field->ts_enabled)
2603           {
2604             CFile* tsFile = CFile::create();
2605             tsFile->duplicateAttributes(file);
2606
2607             // Add variables originating from file and targeted to timeserie file
2608             for (size_t k = 0; k < fileVars.size(); k++)
2609               tsFile->getVirtualVariableGroup()->addChild(fileVars[k]);
2610
2611           
2612             tsFile->name = tsPrefix + "_";
2613             if (!field->name.isEmpty())
2614               tsFile->name.get() += field->name;
2615             else if (field->hasDirectFieldReference()) // We cannot use getBaseFieldReference() just yet
2616               tsFile->name.get() += field->field_ref;
2617             else
2618               tsFile->name.get() += field->getId();
2619
2620             if (!field->ts_split_freq.isEmpty())
2621               tsFile->split_freq = field->ts_split_freq;
2622
2623             CField* tsField = tsFile->addField();
2624             tsField->field_ref = field->getId();
2625
2626             // Add variables originating from file and targeted to timeserie field
2627             for (size_t k = 0; k < fieldVars.size(); k++)
2628               tsField->getVirtualVariableGroup()->addChild(fieldVars[k]);
2629
2630             vars = field->getAllVariables();
2631             for (size_t k = 0; k < vars.size(); k++)
2632             {
2633               CVariable* var = vars[k];
2634
2635               // Add variables originating from field and targeted to timeserie field
2636               if (var->ts_target.isEmpty()
2637                    || var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both)
2638                 tsField->getVirtualVariableGroup()->addChild(var);
2639
2640               // Add variables originating from field and targeted to timeserie file
2641               if (!var->ts_target.isEmpty()
2642                    && (var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both))
2643                 tsFile->getVirtualVariableGroup()->addChild(var);
2644             }
2645
2646             tsFile->solveFieldRefInheritance(true);
2647
2648             if (file->timeseries == CFile::timeseries_attr::exclusive)
2649               field->enabled = false;
2650           }
2651         }
2652
2653         // Finally disable the original file is need be
2654         if (file->timeseries == CFile::timeseries_attr::only)
2655          file->enabled = false;
2656       }
2657     }
2658   }
2659   CATCH_DUMP_ATTR
2660
2661   //! Client side: Send information of reference grid of active fields
2662   void CContext::sendRefGrid(const std::vector<CFile*>& activeFiles)
2663   TRY
2664   {
2665     std::set<pair<StdString,CContextClient*>> gridIds;
2666
2667     int sizeFile = activeFiles.size();
2668     CFile* filePtr(NULL);
2669
2670     // Firstly, find all reference grids of all active fields
2671     for (int i = 0; i < sizeFile; ++i)
2672     {
2673       filePtr = activeFiles[i];
2674       std::vector<CField*> enabledFields = filePtr->getEnabledFields();
2675       int sizeField = enabledFields.size();
2676       for (int numField = 0; numField < sizeField; ++numField)
2677       {
2678         if (0 != enabledFields[numField]->getRelGrid())
2679           gridIds.insert(make_pair(CGrid::get(enabledFields[numField]->getRelGrid())->getId(),enabledFields[numField]->getContextClient()));
2680       }
2681     }
2682
2683     // Create all reference grids on server side
2684     StdString gridDefRoot("grid_definition");
2685     CGridGroup* gridPtr = CGridGroup::get(gridDefRoot);
2686     for (auto it = gridIds.begin(); it != gridIds.end(); ++it)
2687     {
2688       gridPtr->sendCreateChild(it->first,it->second);
2689       CGrid::get(it->first)->sendAllAttributesToServer(it->second);
2690       CGrid::get(it->first)->sendAllDomains(it->second);
2691       CGrid::get(it->first)->sendAllAxis(it->second);
2692       CGrid::get(it->first)->sendAllScalars(it->second);
2693     }
2694   }
2695   CATCH_DUMP_ATTR
2696
2697   //! Client side: Send information of reference domain, axis and scalar of active fields
2698   void CContext::sendRefDomainsAxisScalars(const std::vector<CFile*>& activeFiles)
2699   TRY
2700   {
2701     std::set<pair<StdString,CContextClient*>> domainIds, axisIds, scalarIds;
2702
2703     // Find all reference domain and axis of all active fields
2704     int numEnabledFiles = activeFiles.size();
2705     for (int i = 0; i < numEnabledFiles; ++i)
2706     {
2707       std::vector<CField*> enabledFields = activeFiles[i]->getEnabledFields();
2708       int numEnabledFields = enabledFields.size();
2709       for (int j = 0; j < numEnabledFields; ++j)
2710       {
2711         CContextClient* contextClient=enabledFields[j]->getContextClient() ;
2712         const std::vector<StdString>& prDomAxisScalarId = enabledFields[j]->getRefDomainAxisIds();
2713         if ("" != prDomAxisScalarId[0]) domainIds.insert(make_pair(prDomAxisScalarId[0],contextClient));
2714         if ("" != prDomAxisScalarId[1]) axisIds.insert(make_pair(prDomAxisScalarId[1],contextClient));
2715         if ("" != prDomAxisScalarId[2]) scalarIds.insert(make_pair(prDomAxisScalarId[2],contextClient));
2716       }
2717     }
2718
2719     // Create all reference axis on server side
2720     std::set<StdString>::iterator itDom, itAxis, itScalar;
2721     std::set<StdString>::const_iterator itE;
2722
2723     StdString scalarDefRoot("scalar_definition");
2724     CScalarGroup* scalarPtr = CScalarGroup::get(scalarDefRoot);
2725     
2726     for (auto itScalar = scalarIds.begin(); itScalar != scalarIds.end(); ++itScalar)
2727     {
2728       if (!itScalar->first.empty())
2729       {
2730         scalarPtr->sendCreateChild(itScalar->first,itScalar->second);
2731         CScalar::get(itScalar->first)->sendAllAttributesToServer(itScalar->second);
2732       }
2733     }
2734
2735     StdString axiDefRoot("axis_definition");
2736     CAxisGroup* axisPtr = CAxisGroup::get(axiDefRoot);
2737     
2738     for (auto itAxis = axisIds.begin(); itAxis != axisIds.end(); ++itAxis)
2739     {
2740       if (!itAxis->first.empty())
2741       {
2742         axisPtr->sendCreateChild(itAxis->first, itAxis->second);
2743         CAxis::get(itAxis->first)->sendAllAttributesToServer(itAxis->second);
2744       }
2745     }
2746
2747     // Create all reference domains on server side
2748     StdString domDefRoot("domain_definition");
2749     CDomainGroup* domPtr = CDomainGroup::get(domDefRoot);
2750     
2751     for (auto itDom = domainIds.begin(); itDom != domainIds.end(); ++itDom)
2752     {
2753       if (!itDom->first.empty()) {
2754          domPtr->sendCreateChild(itDom->first, itDom->second);
2755          CDomain::get(itDom->first)->sendAllAttributesToServer(itDom->second);
2756       }
2757     }
2758   }
2759   CATCH_DUMP_ATTR
2760
2761   void CContext::triggerLateFields(void)
2762   TRY
2763   {
2764    for(auto& field : fileInFields_) field->triggerLateField() ;
2765    for(auto& field : couplerInFields_) field->triggerLateField() ;
2766   }
2767   CATCH_DUMP_ATTR
2768
2769   //! Update calendar in each time step
2770   void CContext::updateCalendar(int step)
2771   TRY
2772   {
2773      int prevStep = calendar->getStep();
2774
2775      if (prevStep < step)
2776      {
2777        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2778        {
2779          triggerLateFields();
2780        }
2781
2782        info(50) << "updateCalendar : before : " << calendar->getCurrentDate() << endl;
2783        calendar->update(step);
2784        info(50) << "updateCalendar : after : " << calendar->getCurrentDate() << endl;
2785  #ifdef XIOS_MEMTRACK_LIGHT
2786        info(50) << " Current memory used by XIOS : "<<  MemTrack::getCurrentMemorySize()*1.0/(1024*1024)<<" Mbyte, at timestep "<<step<<" of context "<<this->getId()<<endl ;
2787  #endif
2788
2789        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2790        {
2791          doPostTimestepOperationsForEnabledReadModeFiles();
2792          garbageCollector.invalidate(calendar->getCurrentDate());
2793        }
2794      }
2795      else if (prevStep == step)
2796        info(50) << "updateCalendar: already at step " << step << ", no operation done." << endl;
2797      else // if (prevStep > step)
2798        ERROR("void CContext::updateCalendar(int step)",
2799              << "Illegal calendar update: previous step was " << prevStep << ", new step " << step << "is in the past!")
2800   }
2801   CATCH_DUMP_ATTR
2802
2803   void CContext::initReadFiles(void)
2804   TRY
2805   {
2806      vector<CFile*>::const_iterator it;
2807
2808      for (it=enabledReadModeFiles.begin(); it != enabledReadModeFiles.end(); it++)
2809      {
2810         (*it)->initRead();
2811      }
2812   }
2813   CATCH_DUMP_ATTR
2814
2815   //! Server side: Create header of netcdf file
2816   void CContext::createFileHeader(void)
2817   TRY
2818   {
2819      vector<CFile*>::const_iterator it;
2820
2821      //for (it=enabledFiles.begin(); it != enabledFiles.end(); it++)
2822      for (it=enabledWriteModeFiles.begin(); it != enabledWriteModeFiles.end(); it++)
2823      {
2824         (*it)->initWrite();
2825      }
2826   }
2827   CATCH_DUMP_ATTR
2828
2829   //! Get current context
2830   CContext* CContext::getCurrent(void)
2831   TRY
2832   {
2833     return CObjectFactory::GetObject<CContext>(CObjectFactory::GetCurrentContextId()).get();
2834   }
2835   CATCH
2836
2837   /*!
2838   \brief Set context with an id be the current context
2839   \param [in] id identity of context to be set to current
2840   */
2841   void CContext::setCurrent(const string& id)
2842   TRY
2843   {
2844     CObjectFactory::SetCurrentContextId(id);
2845     CGroupFactory::SetCurrentContextId(id);
2846   }
2847   CATCH
2848
2849  /*!
2850  \brief Create a context with specific id
2851  \param [in] id identity of new context
2852  \return pointer to the new context or already-existed one with identity id
2853  */
2854  CContext* CContext::create(const StdString& id)
2855  TRY
2856  {
2857    CContext::setCurrent(id);
2858
2859    bool hasctxt = CContext::has(id);
2860    CContext* context = CObjectFactory::CreateObject<CContext>(id).get();
2861    getRoot();
2862    if (!hasctxt) CGroupFactory::AddChild(root, context->getShared());
2863
2864#define DECLARE_NODE(Name_, name_) \
2865    C##Name_##Definition::create(C##Name_##Definition::GetDefName());
2866#define DECLARE_NODE_PAR(Name_, name_)
2867#include "node_type.conf"
2868
2869    return (context);
2870  }
2871  CATCH
2872
2873     //! Server side: Receive a message to do some post processing
2874  void CContext::recvRegistry(CEventServer& event)
2875  TRY
2876  {
2877    CBufferIn* buffer=event.subEvents.begin()->buffer;
2878    getCurrent()->recvRegistry(*buffer);
2879  }
2880  CATCH
2881
2882  void CContext::recvRegistry(CBufferIn& buffer)
2883  TRY
2884  {
2885    if (server->intraCommRank==0)
2886    {
2887      CRegistry registry(server->intraComm) ;
2888      registry.fromBuffer(buffer) ;
2889      registryOut->mergeRegistry(registry) ;
2890    }
2891  }
2892  CATCH_DUMP_ATTR
2893
2894  void CContext::sendRegistry(void)
2895  TRY
2896  {
2897    registryOut->hierarchicalGatherRegistry() ;
2898
2899    int nbSrvPools ;
2900    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2901    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2902    else nbSrvPools = 0 ;
2903    CContextClient* contextClientTmp ;
2904
2905    for (int i = 0; i < nbSrvPools; ++i)
2906    {
2907      if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2908      else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2909
2910      CEventClient event(CContext::GetType(), CContext::EVENT_ID_SEND_REGISTRY);
2911      if (contextClientTmp->isServerLeader())
2912      {
2913        CMessage msg ;
2914        if (contextClientTmp->clientRank==0) msg<<*registryOut ;
2915        const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2916        for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2917             event.push(*itRank,1,msg);
2918        contextClientTmp->sendEvent(event);
2919      }
2920      else contextClientTmp->sendEvent(event);
2921    }
2922  }
2923  CATCH_DUMP_ATTR
2924
2925 
2926  void CContext::sendFinalizeClient(CContextClient* contextClient, const string& contextClientId)
2927  TRY
2928  {
2929    CEventClient event(getType(),EVENT_ID_CONTEXT_FINALIZE_CLIENT);
2930    if (contextClient->isServerLeader())
2931    {
2932      CMessage msg;
2933      msg<<contextClientId ;
2934      const std::list<int>& ranks = contextClient->getRanksServerLeader();
2935      for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2936           event.push(*itRank,1,msg);
2937      contextClient->sendEvent(event);
2938    }
2939    else contextClient->sendEvent(event);
2940  }
2941  CATCH_DUMP_ATTR
2942
2943 
2944  void CContext::recvFinalizeClient(CEventServer& event)
2945  TRY
2946  {
2947    CBufferIn* buffer=event.subEvents.begin()->buffer;
2948    string id;
2949    *buffer>>id;
2950    get(id)->recvFinalizeClient(*buffer);
2951  }
2952  CATCH
2953
2954  void CContext::recvFinalizeClient(CBufferIn& buffer)
2955  TRY
2956  {
2957    countChildContextFinalized_++ ;
2958  }
2959  CATCH_DUMP_ATTR
2960
2961
2962
2963
2964 //! Client side: Send a message  announcing that context can receive grid definition from coupling
2965   void CContext::sendCouplerInReady(CContextClient* client)
2966   TRY
2967   {
2968      if (sendCouplerInReady_done_.count(client)!=0) return ;
2969      else sendCouplerInReady_done_.insert(client) ;
2970
2971      CEventClient event(getType(),EVENT_ID_COUPLER_IN_READY);
2972
2973      if (client->isServerLeader())
2974      {
2975        CMessage msg;
2976        msg<<this->getId();
2977        for (auto& rank : client->getRanksServerLeader()) event.push(rank,1,msg);
2978        client->sendEvent(event);
2979      }
2980      else client->sendEvent(event);
2981   }
2982   CATCH_DUMP_ATTR
2983
2984   //! Server side: Receive a message announcing that context can send grid definition for context coupling
2985   void CContext::recvCouplerInReady(CEventServer& event)
2986   TRY
2987   {
2988      CBufferIn* buffer=event.subEvents.begin()->buffer;
2989      getCurrent()->recvCouplerInReady(*buffer);
2990   }
2991   CATCH
2992
2993   //! Server side: Receive a message announcing that context can send grid definition for context coupling
2994   void CContext::recvCouplerInReady(CBufferIn& buffer)
2995   TRY
2996   {
2997      string contextId ;
2998      buffer>>contextId;
2999      couplerInReady_.insert(getCouplerOutClient(contextId)) ;
3000   }
3001   CATCH_DUMP_ATTR
3002
3003
3004
3005
3006
3007 //! Client side: Send a message  announcing that a coupling context have done it closeDefinition, so data can be sent now.
3008   void CContext::sendCouplerInCloseDefinition(CContextClient* client)
3009   TRY
3010   {
3011      if (sendCouplerInCloseDefinition_done_.count(client)!=0) return ;
3012      else sendCouplerInCloseDefinition_done_.insert(client) ;
3013
3014      CEventClient event(getType(),EVENT_ID_COUPLER_IN_CLOSE_DEFINITION);
3015
3016      if (client->isServerLeader())
3017      {
3018        CMessage msg;
3019        msg<<this->getId();
3020        for (auto& rank : client->getRanksServerLeader()) event.push(rank,1,msg);
3021        client->sendEvent(event);
3022      }
3023      else client->sendEvent(event);
3024   }
3025   CATCH_DUMP_ATTR
3026
3027   //! Server side: Receive a message announcing that a coupling context have done it closeDefinition, so data can be sent now.
3028   void CContext::recvCouplerInCloseDefinition(CEventServer& event)
3029   TRY
3030   {
3031      CBufferIn* buffer=event.subEvents.begin()->buffer;
3032      getCurrent()->recvCouplerInCloseDefinition(*buffer);
3033   }
3034   CATCH
3035
3036   //! Server side: Receive a message announcing that a coupling context have done it closeDefinition, so data can be sent now.
3037   void CContext::recvCouplerInCloseDefinition(CBufferIn& buffer)
3038   TRY
3039   {
3040      string contextId ;
3041      buffer>>contextId;
3042      couplerInCloseDefinition_.insert(getCouplerOutClient(contextId)) ;
3043   }
3044   CATCH_DUMP_ATTR
3045
3046
3047
3048
3049//! Client side: Send a message  announcing that a coupling context have done it contextFinalize, so it can also close it own context.
3050   void CContext::sendCouplerInContextFinalized(CContextClient* client)
3051   TRY
3052   {
3053      if (sendCouplerInContextFinalized_done_.count(client)!=0) return ;
3054      else sendCouplerInContextFinalized_done_.insert(client) ;
3055
3056      CEventClient event(getType(),EVENT_ID_COUPLER_IN_CONTEXT_FINALIZED);
3057
3058      if (client->isServerLeader())
3059      {
3060        CMessage msg;
3061        msg<<this->getId();
3062        for (auto& rank : client->getRanksServerLeader()) event.push(rank,1,msg);
3063        client->sendEvent(event);
3064      }
3065      else client->sendEvent(event);
3066   }
3067   CATCH_DUMP_ATTR
3068
3069   //! Server side: Receive a message announcing that a coupling context have done it contextFinalize, so it can also close it own context.
3070   void CContext::recvCouplerInContextFinalized(CEventServer& event)
3071   TRY
3072   {
3073      CBufferIn* buffer=event.subEvents.begin()->buffer;
3074      getCurrent()->recvCouplerInContextFinalized(*buffer);
3075   }
3076   CATCH
3077
3078   //! Server side: Receive a message announcing that a coupling context have done it contextFinalize, so it can also close it own context.
3079   void CContext::recvCouplerInContextFinalized(CBufferIn& buffer)
3080   TRY
3081   {
3082      string contextId ;
3083      buffer>>contextId;
3084      couplerInContextFinalized_.insert(getCouplerOutClient(contextId)) ;
3085   }
3086   CATCH_DUMP_ATTR
3087
3088
3089
3090
3091  /*!
3092  * \fn bool CContext::isFinalized(void)
3093  * Context is finalized if it received context post finalize event.
3094  */
3095  bool CContext::isFinalized(void)
3096  TRY
3097  {
3098    return finalized;
3099  }
3100  CATCH_DUMP_ATTR
3101  ///--------------------------------------------------------------
3102  StdString CContext::dumpClassAttributes(void)
3103  {
3104    StdString str;
3105    str.append("enabled files=\"");
3106    int size = this->enabledFiles.size();
3107    for (int i = 0; i < size; ++i)
3108    {
3109      str.append(enabledFiles[i]->getId());
3110      str.append(" ");
3111    }
3112    str.append("\"");
3113    return str;
3114  }
3115
3116} // namespace xios
Note: See TracBrowser for help on using the repository browser.