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

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

Some update on XIOS_COUPLING branch...

YM

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