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

Last change on this file since 1871 was 1871, 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.6 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   
1094    // workflow endpoint => sent to IO/SERVER
1095    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1096    {
1097      for(auto field : fileOutField) 
1098      {
1099        field->connectToFileServer(garbageCollector) ; // connect the field to server filter
1100        field->computeGridIndexToFileServer() ; // compute grid index for transfer to the server context
1101      }
1102      setClientServerBuffer(fileOutField, true) ; // set context
1103      for(auto field : fileOutField) field->sendFieldToFileServer() ;
1104    }
1105
1106    // workflow endpoint => write to file
1107    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)
1108    {
1109      for(auto field : fileOutField) 
1110      {
1111        field->connectToFileWriter(garbageCollector) ; // connect the field to server filter
1112      }
1113    }
1114   
1115    // workflow endpoint => Send data from server to client
1116    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER)
1117    {
1118      // no filter to send data from server to client => to be implemented (reading case)
1119    }
1120
1121    // workflow endpoint => sent to model on client side
1122    if (serviceType_==CServicesManager::CLIENT)
1123    {
1124      for(auto field : fieldWithReadAccess) field->connectToModelOutput(garbageCollector) ;
1125    }
1126
1127
1128    // workflow startpoint => data from model
1129    if (serviceType_==CServicesManager::CLIENT)
1130    {
1131      for(auto field : fieldModelIn) 
1132      {
1133        field->connectToModelInput(garbageCollector) ; // connect the field to server filter
1134        // grid index will be computed on the fly
1135      }
1136    }
1137   
1138    // workflow startpoint => data from client on server side
1139    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER || serviceType_==CServicesManager::OUT_SERVER)
1140    {
1141      for(auto field : fieldModelIn) 
1142      {
1143        field->connectToClientInput(garbageCollector) ; // connect the field to server filter
1144      }
1145    }
1146
1147     // workflow startpoint => data from server on client side
1148    if (serviceType_==CServicesManager::CLIENT)
1149    {
1150      for(auto field : fileInField) 
1151      {
1152        field->connectToServerInput(garbageCollector) ; // connect the field to server filter
1153        field->computeGridIndexToFileServer() ; // compute grid index for transfer to the server context
1154        field->sendFieldToFileServer() ;
1155      }
1156    }
1157
1158    // workflow startpoint => data read from file on server side
1159    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::GATHERER)
1160    {
1161      // no filter for reading data from file => to be implemented
1162    }
1163
1164
1165    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1166    if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER)  createFileHeader();
1167    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1168   
1169
1170
1171
1172    return ;
1173    // For now, only read files with client and only one level server
1174    // if (hasClient && !hasServer) findEnabledReadModeFiles();     
1175
1176    // Find all enabled fields of each file     
1177    findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
1178    findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
1179
1180    // For now, only read files with client and only one level server
1181    // if (hasClient && !hasServer)
1182    //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
1183
1184    if (serviceType_==CServicesManager::CLIENT)
1185    {
1186      initReadFiles();
1187      // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
1188      this->readAttributesOfEnabledFieldsInReadModeFiles();
1189    }
1190
1191    // Only search and rebuild all reference objects of enable fields, don't transform
1192    this->solveOnlyRefOfEnabledFields();
1193
1194    // Search and rebuild all reference object of enabled fields, and transform
1195    this->solveAllRefOfEnabledFieldsAndTransform();
1196
1197    // Find all fields with read access from the public API
1198    if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
1199    // and solve the all reference for them
1200    if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
1201
1202    isPostProcessed = true;
1203
1204
1205
1206    // Distribute files between secondary servers according to the data size
1207    distributeFiles(this->enabledWriteModeFiles);
1208
1209    // Check grid and calculate its distribution
1210    checkGridEnabledFields();
1211
1212    setClientServerBuffer(client, (serviceType_==CServicesManager::CLIENT) ) ;
1213    for (int i = 0; i < clientPrimServer.size(); ++i)
1214         setClientServerBuffer(clientPrimServer[i], true);
1215
1216   
1217    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)
1218    { 
1219      if (serviceType_==CServicesManager::GATHERER)
1220      { 
1221        for (auto it=clientPrimServer.begin(); it!=clientPrimServer.end();++it) 
1222        {
1223          this->sendAllAttributesToServer(*it); // Send all attributes of current context to server
1224          CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(*it); // Send all attributes of current calendar
1225        }
1226      }
1227      else 
1228      {
1229        this->sendAllAttributesToServer(client);   // Send all attributes of current context to server
1230        CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current calendar
1231      }
1232
1233      // We have enough information to send to server
1234      // First of all, send all enabled files
1235      sendEnabledFiles(this->enabledWriteModeFiles);
1236      // We only use server-level 1 (for now) to read data
1237      if (serviceType_==CServicesManager::CLIENT)  sendEnabledFiles(this->enabledReadModeFiles);
1238
1239      // Then, send all enabled fields     
1240      sendEnabledFieldsInFiles(this->enabledWriteModeFiles);
1241     
1242      if (serviceType_==CServicesManager::CLIENT) sendEnabledFieldsInFiles(this->enabledReadModeFiles);
1243
1244      // Then, check whether we have domain_ref, axis_ref or scalar_ref attached to the enabled fields
1245      // If any, so send them to server
1246      sendRefDomainsAxisScalars(this->enabledWriteModeFiles); 
1247     
1248      if (serviceType_==CServicesManager::CLIENT) sendRefDomainsAxisScalars(this->enabledReadModeFiles);       
1249
1250      // Check whether enabled fields have grid_ref, if any, send this info to server
1251      sendRefGrid(this->enabledFiles);
1252      // This code may be useful in the future when we want to seperate completely read and write
1253      // sendRefGrid(this->enabledWriteModeFiles);
1254      // if (!hasServer)
1255      //   sendRefGrid(this->enabledReadModeFiles);
1256     
1257      // A grid of enabled fields composed of several components which must be checked then their
1258      // checked attributes should be sent to server
1259      sendGridComponentEnabledFieldsInFiles(this->enabledFiles); // This code can be seperated in two (one for reading, another for writing)
1260
1261      // We have a xml tree on the server side and now, it should be also processed
1262      sendPostProcessing();
1263       
1264      // Finally, we send information of grid itself to server
1265      sendGridEnabledFieldsInFiles(this->enabledWriteModeFiles);       
1266     
1267      if (serviceType_==CServicesManager::CLIENT) sendGridEnabledFieldsInFiles(this->enabledReadModeFiles);       
1268    }
1269    allProcessed = true;
1270
1271
1272    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1273
1274    // There are some processings that should be done after all of above. For example: check mask or index
1275    this->buildFilterGraphOfEnabledFields();
1276   
1277     if (serviceType_==CServicesManager::CLIENT)
1278    {
1279      buildFilterGraphOfFieldsWithReadAccess();
1280      postProcessFilterGraph(); // For coupling in, modify this later
1281    }
1282   
1283    checkGridEnabledFields();   
1284
1285    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1286    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1287
1288    // Nettoyage de l'arborescence
1289    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1290
1291    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1292    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1293   
1294    CTimer::get("Context : close definition").suspend() ;
1295  }
1296  CATCH_DUMP_ATTR
1297
1298 /*!
1299  * Send context attribute and calendar to file server, it must be done once by context file server
1300  * \param[in] client : context client to send   
1301  */ 
1302  void CContext::sendContextToFileServer(CContextClient* client)
1303  {
1304    if (sendToFileServer_done_.count(client)!=0) return ;
1305    else sendToFileServer_done_.insert(client) ;
1306   
1307    this->sendAllAttributesToServer(client); // Send all attributes of current context to server
1308    CCalendarWrapper::get(CCalendarWrapper::GetDefName())->sendAllAttributesToServer(client); // Send all attributes of current cale
1309  }
1310
1311  // ym obsolete now to be removed
1312   void CContext::closeDefinition_old(void)
1313   TRY
1314   {
1315     CTimer::get("Context : close definition").resume() ;
1316   
1317    //
1318    postProcessingGlobalAttributes();
1319
1320    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendPostProcessingGlobalAttributes();
1321
1322    // There are some processings that should be done after all of above. For example: check mask or index
1323    this->buildFilterGraphOfEnabledFields();
1324   
1325     if (serviceType_==CServicesManager::CLIENT)
1326    {
1327      buildFilterGraphOfFieldsWithReadAccess();
1328      postProcessFilterGraph(); // For coupling in, modify this later
1329    }
1330   
1331    checkGridEnabledFields();   
1332
1333    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendProcessingGridOfEnabledFields();
1334    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) this->sendCloseDefinition();
1335
1336    // Nettoyage de l'arborescence
1337    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) CleanTree(); // Only on client side??
1338
1339    if (serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER) sendCreateFileHeader();
1340    if (serviceType_==CServicesManager::CLIENT) startPrefetchingOfEnabledReadModeFiles();
1341   
1342    CTimer::get("Context : close definition").suspend() ;
1343   }
1344   CATCH_DUMP_ATTR
1345
1346   vector<CField*> CContext::findAllEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1347   TRY
1348   {
1349     vector<CField*> fields ;
1350     for (unsigned int i = 0; i < activeFiles.size(); i++)
1351     {
1352        const vector<CField*>&& field=activeFiles[i]->getEnabledFields() ;
1353        fields.insert(fields.end(),field.begin(),field.end());
1354     }
1355     return fields ;
1356   }
1357   CATCH_DUMP_ATTR
1358
1359   vector<CField*> CContext::findAllEnabledFieldsInFileOut(const std::vector<CFile*>& activeFiles)
1360   TRY
1361   {
1362     vector<CField*> fields ;
1363     for(auto file : activeFiles)
1364     {
1365        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1366        for(auto field : fieldList) field->setFileOut(file) ;
1367        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1368     }
1369     return fields ;
1370   }
1371   CATCH_DUMP_ATTR
1372
1373   vector<CField*> CContext::findAllEnabledFieldsInFileIn(const std::vector<CFile*>& activeFiles)
1374   TRY
1375   {
1376     vector<CField*> fields ;
1377     for(auto file : activeFiles)
1378     {
1379        const vector<CField*>&& fieldList=file->getEnabledFields() ;
1380        for(auto field : fieldList) field->setFileIn(file) ;
1381        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1382     }
1383     return fields ;
1384   }
1385   CATCH_DUMP_ATTR
1386
1387   vector<CField*> CContext::findAllEnabledFieldsCouplerOut(const std::vector<CCouplerOut*>& activeCouplerOut)
1388   TRY
1389   {
1390     vector<CField*> fields ;
1391     for (auto couplerOut :activeCouplerOut)
1392     {
1393        const vector<CField*>&& fieldList=couplerOut->getEnabledFields() ;
1394        for(auto field : fieldList) field->setCouplerOut(couplerOut) ;
1395        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1396     }
1397     return fields ;
1398   }
1399   CATCH_DUMP_ATTR
1400
1401   vector<CField*> CContext::findAllEnabledFieldsCouplerIn(const std::vector<CCouplerIn*>& activeCouplerIn)
1402   TRY
1403   {
1404     vector<CField*> fields ;
1405     for (auto couplerIn :activeCouplerIn)
1406     {
1407        const vector<CField*>&& fieldList=couplerIn->getEnabledFields() ;
1408        for(auto field : fieldList) field->setCouplerIn(couplerIn) ;
1409        fields.insert(fields.end(),fieldList.begin(),fieldList.end());
1410     }
1411     return fields ;
1412   }
1413   CATCH_DUMP_ATTR
1414
1415
1416
1417   void CContext::readAttributesOfEnabledFieldsInReadModeFiles()
1418   TRY
1419   {
1420      for (unsigned int i = 0; i < this->enabledReadModeFiles.size(); ++i)
1421        (void)this->enabledReadModeFiles[i]->readAttributesOfEnabledFieldsInReadMode();
1422   }
1423   CATCH_DUMP_ATTR
1424
1425   void CContext::sendGridComponentEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1426   TRY
1427   {
1428     int size = activeFiles.size();
1429     for (int i = 0; i < size; ++i)
1430     {       
1431       activeFiles[i]->sendGridComponentOfEnabledFields();
1432     }
1433   }
1434   CATCH_DUMP_ATTR
1435
1436   /*!
1437      Send active (enabled) fields in file from a client to others
1438      \param [in] activeFiles files contains enabled fields to send
1439   */
1440   void CContext::sendGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1441   TRY
1442   {
1443     int size = activeFiles.size();
1444     for (int i = 0; i < size; ++i)
1445     {       
1446       activeFiles[i]->sendGridOfEnabledFields();
1447     }
1448   }
1449   CATCH_DUMP_ATTR
1450
1451   void CContext::checkGridEnabledFields()
1452   TRY
1453   {
1454     int size = enabledFiles.size();
1455     for (int i = 0; i < size; ++i)
1456     {
1457       enabledFiles[i]->checkGridOfEnabledFields();       
1458     }
1459
1460     size = enabledCouplerOut.size();
1461     for (int i = 0; i < size; ++i)
1462     {
1463       enabledCouplerOut[i]->checkGridOfEnabledFields();       
1464     }
1465   }
1466   CATCH_DUMP_ATTR
1467
1468   /*!
1469      Check grid of active (enabled) fields in file
1470      \param [in] activeFiles files contains enabled fields whose grid needs checking
1471   */
1472   void CContext::checkGridEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
1473   TRY
1474   {
1475     int size = activeFiles.size();
1476     for (int i = 0; i < size; ++i)
1477     {
1478       activeFiles[i]->checkGridOfEnabledFields();       
1479     }
1480   }
1481   CATCH_DUMP_ATTR
1482
1483    /*!
1484      Go up the hierachical tree via field_ref and do check of attributes of fields
1485      This can be done in a client then all computed information will be sent from this client to others
1486      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1487   */
1488   void CContext::solveOnlyRefOfEnabledFields(void)
1489   TRY
1490   {
1491     int size = this->enabledFiles.size();
1492     for (int i = 0; i < size; ++i)
1493     {
1494       this->enabledFiles[i]->solveOnlyRefOfEnabledFields();
1495     }
1496
1497     for (int i = 0; i < size; ++i)
1498     {
1499       this->enabledFiles[i]->generateNewTransformationGridDest();
1500     }
1501
1502     size = this->enabledCouplerOut.size();
1503     for (int i = 0; i < size; ++i)
1504     {
1505       this->enabledCouplerOut[i]->solveOnlyRefOfEnabledFields();
1506     }
1507
1508     for (int i = 0; i < size; ++i)
1509     {
1510       this->enabledCouplerOut[i]->generateNewTransformationGridDest();
1511     }
1512   }
1513   CATCH_DUMP_ATTR
1514
1515    /*!
1516      Go up the hierachical tree via field_ref and do check of attributes of fields.
1517      The transformation can be done in this step.
1518      All computed information will be sent from this client to others.
1519      \param [in] sendToServer Flag to indicate whether calculated information will be sent
1520   */
1521   void CContext::solveAllRefOfEnabledFieldsAndTransform(void)
1522   TRY
1523   {
1524     int size = this->enabledFiles.size();
1525     for (int i = 0; i < size; ++i)
1526     {
1527       this->enabledFiles[i]->solveAllRefOfEnabledFieldsAndTransform();
1528     }
1529
1530     size = this->enabledCouplerOut.size();
1531     for (int i = 0; i < size; ++i)
1532     {
1533       this->enabledCouplerOut[i]->solveAllRefOfEnabledFieldsAndTransform();
1534     }
1535
1536   }
1537   CATCH_DUMP_ATTR
1538
1539   void CContext::buildFilterGraphOfEnabledFields()
1540   TRY
1541   {
1542     int size = this->enabledFiles.size();
1543     for (int i = 0; i < size; ++i)
1544     {
1545       this->enabledFiles[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1546     }
1547
1548     size = this->enabledCouplerOut.size();
1549     for (int i = 0; i < size; ++i)
1550     {
1551       this->enabledCouplerOut[i]->buildFilterGraphOfEnabledFields(garbageCollector);
1552     }
1553   }
1554   CATCH_DUMP_ATTR
1555
1556   void CContext::postProcessFilterGraph()
1557   TRY
1558   {
1559     int size = enabledFiles.size();
1560     for (int i = 0; i < size; ++i)
1561     {
1562        enabledFiles[i]->postProcessFilterGraph();
1563     }
1564   }
1565   CATCH_DUMP_ATTR
1566
1567   void CContext::startPrefetchingOfEnabledReadModeFiles()
1568   TRY
1569   {
1570     int size = enabledReadModeFiles.size();
1571     for (int i = 0; i < size; ++i)
1572     {
1573        enabledReadModeFiles[i]->prefetchEnabledReadModeFields();
1574     }
1575   }
1576   CATCH_DUMP_ATTR
1577
1578   void CContext::doPreTimestepOperationsForEnabledReadModeFiles()
1579   TRY
1580   {
1581     int size = enabledReadModeFiles.size();
1582     for (int i = 0; i < size; ++i)
1583     {
1584        enabledReadModeFiles[i]->doPreTimestepOperationsForEnabledReadModeFields();
1585     }
1586   }
1587   CATCH_DUMP_ATTR
1588
1589   void CContext::doPostTimestepOperationsForEnabledReadModeFiles()
1590   TRY
1591   {
1592     int size = enabledReadModeFiles.size();
1593     for (int i = 0; i < size; ++i)
1594     {
1595        enabledReadModeFiles[i]->doPostTimestepOperationsForEnabledReadModeFields();
1596     }
1597   }
1598   CATCH_DUMP_ATTR
1599
1600  void CContext::findFieldsWithReadAccess(void)
1601  TRY
1602  {
1603    fieldsWithReadAccess_.clear();
1604    const vector<CField*> allFields = CField::getAll();
1605    for (size_t i = 0; i < allFields.size(); ++i)
1606    {
1607      CField* field = allFields[i];
1608      if (!field->read_access.isEmpty() && field->read_access && (field->enabled.isEmpty() || field->enabled))
1609      {
1610        fieldsWithReadAccess_.push_back(field);
1611        field->setModelOut() ;
1612      }
1613    }
1614  }
1615  CATCH_DUMP_ATTR
1616
1617  void CContext::solveAllRefOfFieldsWithReadAccess()
1618  TRY
1619  {
1620    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1621      fieldsWithReadAccess_[i]->solveAllReferenceEnabledField(false);
1622  }
1623  CATCH_DUMP_ATTR
1624
1625  void CContext::buildFilterGraphOfFieldsWithReadAccess()
1626  TRY
1627  {
1628    for (size_t i = 0; i < fieldsWithReadAccess_.size(); ++i)
1629      fieldsWithReadAccess_[i]->buildFilterGraph(garbageCollector, true);
1630  }
1631  CATCH_DUMP_ATTR
1632
1633   void CContext::solveAllInheritance(bool apply)
1634   TRY
1635   {
1636     // Résolution des héritages descendants (càd des héritages de groupes)
1637     // pour chacun des contextes.
1638      solveDescInheritance(apply);
1639
1640     // Résolution des héritages par référence au niveau des fichiers.
1641      const vector<CFile*> allFiles=CFile::getAll();
1642      const vector<CCouplerIn*> allCouplerIn=CCouplerIn::getAll();
1643      const vector<CCouplerOut*> allCouplerOut=CCouplerOut::getAll();
1644      const vector<CGrid*> allGrids= CGrid::getAll();
1645
1646      if (serviceType_==CServicesManager::CLIENT)
1647      {
1648        for (unsigned int i = 0; i < allFiles.size(); i++)
1649          allFiles[i]->solveFieldRefInheritance(apply);
1650
1651        for (unsigned int i = 0; i < allCouplerIn.size(); i++)
1652          allCouplerIn[i]->solveFieldRefInheritance(apply);
1653
1654        for (unsigned int i = 0; i < allCouplerOut.size(); i++)
1655          allCouplerOut[i]->solveFieldRefInheritance(apply);
1656      }
1657
1658      unsigned int vecSize = allGrids.size();
1659      unsigned int i = 0;
1660      for (i = 0; i < vecSize; ++i)
1661        allGrids[i]->solveElementsRefInheritance(apply);
1662
1663   }
1664  CATCH_DUMP_ATTR
1665
1666   void CContext::findEnabledFiles(void)
1667   TRY
1668   {
1669      const std::vector<CFile*> allFiles = CFile::getAll();
1670      const CDate& initDate = calendar->getInitDate();
1671
1672      for (unsigned int i = 0; i < allFiles.size(); i++)
1673         if (!allFiles[i]->enabled.isEmpty()) // Si l'attribut 'enabled' est défini.
1674         {
1675            if (allFiles[i]->enabled.getValue()) // Si l'attribut 'enabled' est fixé à vrai.
1676            {
1677              if (allFiles[i]->output_freq.isEmpty())
1678              {
1679                 ERROR("CContext::findEnabledFiles()",
1680                     << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1681                     <<" \".")
1682              }
1683              if ((initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1684              {
1685                error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1686                    << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1687                    <<"\" is less than the time step. File will not be written."<<endl;
1688              }
1689              else
1690               enabledFiles.push_back(allFiles[i]);
1691            }
1692         }
1693         else
1694         {
1695           if (allFiles[i]->output_freq.isEmpty())
1696           {
1697              ERROR("CContext::findEnabledFiles()",
1698                  << "Mandatory attribute output_freq must be defined for file \""<<allFiles[i]->getFileOutputName()
1699                  <<" \".")
1700           }
1701           if ( (initDate + allFiles[i]->output_freq.getValue()) < (initDate + this->getCalendar()->getTimeStep()))
1702           {
1703             error(0)<<"WARNING: void CContext::findEnabledFiles()"<<endl
1704                 << "Output frequency in file \""<<allFiles[i]->getFileOutputName()
1705                 <<"\" is less than the time step. File will not be written."<<endl;
1706           }
1707           else
1708             enabledFiles.push_back(allFiles[i]); // otherwise true by default
1709         }
1710
1711      if (enabledFiles.size() == 0)
1712         DEBUG(<<"Aucun fichier ne va être sorti dans le contexte nommé \""
1713               << getId() << "\" !");
1714
1715   }
1716   CATCH_DUMP_ATTR
1717
1718   void CContext::findEnabledCouplerIn(void)
1719   TRY
1720   {
1721      const std::vector<CCouplerIn*> allCouplerIn = CCouplerIn::getAll();
1722      bool enabled ;
1723      for (size_t i = 0; i < allCouplerIn.size(); i++)
1724      {
1725        if (allCouplerIn[i]->enabled.isEmpty()) enabled=true ;
1726        else enabled=allCouplerIn[i]->enabled ;
1727        if (enabled) enabledCouplerIn.push_back(allCouplerIn[i]) ;
1728      }
1729   }
1730   CATCH_DUMP_ATTR
1731
1732   void CContext::findEnabledCouplerOut(void)
1733   TRY
1734   {
1735      const std::vector<CCouplerOut*> allCouplerOut = CCouplerOut::getAll();
1736      bool enabled ;
1737      for (size_t i = 0; i < allCouplerOut.size(); i++)
1738      {
1739        if (allCouplerOut[i]->enabled.isEmpty()) enabled=true ;
1740        else enabled=allCouplerOut[i]->enabled ;
1741        if (enabled) enabledCouplerOut.push_back(allCouplerOut[i]) ;
1742      }
1743   }
1744   CATCH_DUMP_ATTR
1745
1746
1747
1748
1749   void CContext::distributeFiles(const vector<CFile*>& files)
1750   TRY
1751   {
1752     bool distFileMemory=false ;
1753     distFileMemory=CXios::getin<bool>("server2_dist_file_memory", distFileMemory);
1754
1755     if (distFileMemory) distributeFileOverMemoryBandwith(files) ;
1756     else distributeFileOverBandwith(files) ;
1757   }
1758   CATCH_DUMP_ATTR
1759
1760   void CContext::distributeFileOverBandwith(const vector<CFile*>& files)
1761   TRY
1762   {
1763     double eps=std::numeric_limits<double>::epsilon()*10 ;
1764     
1765     std::ofstream ofs(("distribute_file_"+getId()+".dat").c_str(), std::ofstream::out);
1766     int nbPools = clientPrimServer.size();
1767
1768     // (1) Find all enabled files in write mode
1769     // for (int i = 0; i < this->enabledFiles.size(); ++i)
1770     // {
1771     //   if (enabledFiles[i]->mode.isEmpty() || (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
1772     //    enabledWriteModeFiles.push_back(enabledFiles[i]);
1773     // }
1774
1775     // (2) Estimate the data volume for each file
1776     int size = files.size();
1777     std::vector<std::pair<double, CFile*> > dataSizeMap;
1778     double dataPerPool = 0;
1779     int nfield=0 ;
1780     ofs<<size<<endl ;
1781     for (size_t i = 0; i < size; ++i)
1782     {
1783       CFile* file = files[i];
1784       ofs<<file->getId()<<endl ;
1785       StdSize dataSize=0;
1786       std::vector<CField*> enabledFields = file->getEnabledFields();
1787       size_t numEnabledFields = enabledFields.size();
1788       ofs<<numEnabledFields<<endl ;
1789       for (size_t j = 0; j < numEnabledFields; ++j)
1790       {
1791         dataSize += enabledFields[j]->getGlobalWrittenSize() ;
1792         ofs<<enabledFields[j]->getGrid()->getId()<<endl ;
1793         ofs<<enabledFields[j]->getGlobalWrittenSize()<<endl ;
1794       }
1795       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
1796       double dataSizeSec= dataSize/ outFreqSec;
1797       ofs<<dataSizeSec<<endl ;
1798       nfield++ ;
1799// add epsilon*nField to dataSizeSec in order to  preserve reproductive ordering when sorting
1800       dataSizeMap.push_back(make_pair(dataSizeSec + dataSizeSec * eps * nfield , file));
1801       dataPerPool += dataSizeSec;
1802     }
1803     dataPerPool /= nbPools;
1804     std::sort(dataSizeMap.begin(), dataSizeMap.end());
1805
1806     // (3) Assign contextClient to each enabled file
1807
1808     std::multimap<double,int> poolDataSize ;
1809// multimap is not garanty to preserve stable sorting in c++98 but it seems it does for c++11
1810
1811     int j;
1812     double dataSize ;
1813     for (j = 0 ; j < nbPools ; ++j) poolDataSize.insert(std::pair<double,int>(0.,j)) ; 
1814             
1815     for (int i = dataSizeMap.size()-1; i >= 0; --i)
1816     {
1817       dataSize=(*poolDataSize.begin()).first ;
1818       j=(*poolDataSize.begin()).second ;
1819       dataSizeMap[i].second->setContextClient(clientPrimServer[j]);
1820       dataSize+=dataSizeMap[i].first;
1821       poolDataSize.erase(poolDataSize.begin()) ;
1822       poolDataSize.insert(std::pair<double,int>(dataSize,j)) ; 
1823     }
1824
1825     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 ;
1826   }
1827   CATCH_DUMP_ATTR
1828
1829   void CContext::distributeFileOverMemoryBandwith(const vector<CFile*>& filesList)
1830   TRY
1831   {
1832     int nbPools = clientPrimServer.size();
1833     double ratio=0.5 ;
1834     ratio=CXios::getin<double>("server2_dist_file_memory_ratio", ratio);
1835
1836     int nFiles = filesList.size();
1837     vector<SDistFile> files(nFiles);
1838     vector<SDistGrid> grids;
1839     map<string,int> gridMap ;
1840     string gridId; 
1841     int gridIndex=0 ;
1842
1843     for (size_t i = 0; i < nFiles; ++i)
1844     {
1845       StdSize dataSize=0;
1846       CFile* file = filesList[i];
1847       std::vector<CField*> enabledFields = file->getEnabledFields();
1848       size_t numEnabledFields = enabledFields.size();
1849
1850       files[i].id_=file->getId() ;
1851       files[i].nbGrids_=numEnabledFields;
1852       files[i].assignedGrid_ = new int[files[i].nbGrids_] ;
1853         
1854       for (size_t j = 0; j < numEnabledFields; ++j)
1855       {
1856         gridId=enabledFields[j]->getGrid()->getId() ;
1857         if (gridMap.find(gridId)==gridMap.end())
1858         {
1859            gridMap[gridId]=gridIndex  ;
1860            SDistGrid newGrid; 
1861            grids.push_back(newGrid) ;
1862            gridIndex++ ;
1863         }
1864         files[i].assignedGrid_[j]=gridMap[gridId] ;
1865         grids[files[i].assignedGrid_[j]].size_=enabledFields[j]->getGlobalWrittenSize() ;
1866         dataSize += enabledFields[j]->getGlobalWrittenSize() ; // usefull
1867       }
1868       double outFreqSec = (Time)(calendar->getCurrentDate()+file->output_freq)-(Time)(calendar->getCurrentDate()) ;
1869       files[i].bandwith_= dataSize/ outFreqSec ;
1870     }
1871
1872     double bandwith=0 ;
1873     double memory=0 ;
1874   
1875     for(int i=0; i<nFiles; i++)  bandwith+=files[i].bandwith_ ;
1876     for(int i=0; i<nFiles; i++)  files[i].bandwith_ = files[i].bandwith_/bandwith * ratio ;
1877
1878     for(int i=0; i<grids.size(); i++)  memory+=grids[i].size_ ;
1879     for(int i=0; i<grids.size(); i++)  grids[i].size_ = grids[i].size_ / memory * (1.0-ratio) ;
1880       
1881     distributeFileOverServer2(nbPools, grids.size(), &grids[0], nFiles, &files[0]) ;
1882
1883     vector<double> memorySize(nbPools,0.) ;
1884     vector< set<int> > serverGrids(nbPools) ;
1885     vector<double> bandwithSize(nbPools,0.) ;
1886       
1887     for (size_t i = 0; i < nFiles; ++i)
1888     {
1889       bandwithSize[files[i].assignedServer_] += files[i].bandwith_* bandwith /ratio ;
1890       for(int j=0 ; j<files[i].nbGrids_;j++)
1891       {
1892         if (serverGrids[files[i].assignedServer_].find(files[i].assignedGrid_[j]) == serverGrids[files[i].assignedServer_].end())
1893         {
1894           memorySize[files[i].assignedServer_]+= grids[files[i].assignedGrid_[j]].size_ * memory / (1.0-ratio);
1895           serverGrids[files[i].assignedServer_].insert(files[i].assignedGrid_[j]) ;
1896         }
1897       }
1898       filesList[i]->setContextClient(clientPrimServer[files[i].assignedServer_]) ;
1899       delete [] files[i].assignedGrid_ ;
1900     }
1901
1902     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 ;
1903     for (int i = 0; i < nbPools; ++i) info(100)<<"Pool server level2 "<<i<<"   assigned grid memory "<<memorySize[i]*100/1024./1024.<<" Mb"<<endl ;
1904
1905   }
1906   CATCH_DUMP_ATTR
1907
1908   /*!
1909      Find all files in write mode
1910   */
1911   void CContext::findEnabledWriteModeFiles(void)
1912   TRY
1913   {
1914     int size = this->enabledFiles.size();
1915     for (int i = 0; i < size; ++i)
1916     {
1917       if (enabledFiles[i]->mode.isEmpty() || 
1918          (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::write ))
1919        enabledWriteModeFiles.push_back(enabledFiles[i]);
1920     }
1921   }
1922   CATCH_DUMP_ATTR
1923
1924   /*!
1925      Find all files in read mode
1926   */
1927   void CContext::findEnabledReadModeFiles(void)
1928   TRY
1929   {
1930     int size = this->enabledFiles.size();
1931     for (int i = 0; i < size; ++i)
1932     {
1933       if (!enabledFiles[i]->mode.isEmpty() && enabledFiles[i]->mode.getValue() == CFile::mode_attr::read)
1934        enabledReadModeFiles.push_back(enabledFiles[i]);
1935     }
1936   }
1937   CATCH_DUMP_ATTR
1938
1939   void CContext::closeAllFile(void)
1940   TRY
1941   {
1942     std::vector<CFile*>::const_iterator
1943            it = this->enabledFiles.begin(), end = this->enabledFiles.end();
1944
1945     for (; it != end; it++)
1946     {
1947       info(30)<<"Closing File : "<<(*it)->getId()<<endl;
1948       (*it)->close();
1949     }
1950   }
1951   CATCH_DUMP_ATTR
1952
1953   /*!
1954   \brief Dispatch event received from client
1955      Whenever a message is received in buffer of server, it will be processed depending on
1956   its event type. A new event type should be added in the switch list to make sure
1957   it processed on server side.
1958   \param [in] event: Received message
1959   */
1960   bool CContext::dispatchEvent(CEventServer& event)
1961   TRY
1962   {
1963
1964      if (SuperClass::dispatchEvent(event)) return true;
1965      else
1966      {
1967        switch(event.type)
1968        {
1969           case EVENT_ID_CLOSE_DEFINITION :
1970             recvCloseDefinition(event);
1971             return true;
1972             break;
1973           case EVENT_ID_UPDATE_CALENDAR:
1974             recvUpdateCalendar(event);
1975             return true;
1976             break;
1977           case EVENT_ID_CREATE_FILE_HEADER :
1978             recvCreateFileHeader(event);
1979             return true;
1980             break;
1981           case EVENT_ID_POST_PROCESS:
1982             recvPostProcessing(event);
1983             return true;
1984            case EVENT_ID_SEND_REGISTRY:
1985             recvRegistry(event);
1986             return true;
1987             break;
1988            case EVENT_ID_POST_PROCESS_GLOBAL_ATTRIBUTES:
1989             recvPostProcessingGlobalAttributes(event);
1990             return true;
1991             break;
1992            case EVENT_ID_PROCESS_GRID_ENABLED_FIELDS:
1993             recvProcessingGridOfEnabledFields(event);
1994             return true;
1995             break;
1996           default :
1997             ERROR("bool CContext::dispatchEvent(CEventServer& event)",
1998                    <<"Unknown Event");
1999           return false;
2000         }
2001      }
2002   }
2003   CATCH
2004
2005   //! Client side: Send a message to server to make it close
2006   void CContext::sendCloseDefinition(void)
2007   TRY
2008   {
2009    int nbSrvPools ;
2010    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2011    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2012    else nbSrvPools = 0 ;
2013    CContextClient* contextClientTmp ;
2014
2015    for (int i = 0; i < nbSrvPools; ++i)
2016     {
2017       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2018       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2019       CEventClient event(getType(),EVENT_ID_CLOSE_DEFINITION);
2020       if (contextClientTmp->isServerLeader())
2021       {
2022         CMessage msg;
2023         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2024         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2025           event.push(*itRank,1,msg);
2026         contextClientTmp->sendEvent(event);
2027       }
2028       else contextClientTmp->sendEvent(event);
2029     }
2030   }
2031   CATCH_DUMP_ATTR
2032
2033   //! Server side: Receive a message of client announcing a context close
2034   void CContext::recvCloseDefinition(CEventServer& event)
2035   TRY
2036   {
2037      CBufferIn* buffer=event.subEvents.begin()->buffer;
2038      getCurrent()->closeDefinition();
2039   }
2040   CATCH
2041
2042   //! Client side: Send a message to update calendar in each time step
2043   void CContext::sendUpdateCalendar(int step)
2044   TRY
2045   {
2046    int nbSrvPools ;
2047    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2048    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2049    else nbSrvPools = 0 ;
2050    CContextClient* contextClientTmp ;
2051
2052    for (int i = 0; i < nbSrvPools; ++i)
2053    {
2054       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2055       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2056       CEventClient event(getType(),EVENT_ID_UPDATE_CALENDAR);
2057
2058         if (contextClientTmp->isServerLeader())
2059         {
2060           CMessage msg;
2061           msg<<step;
2062           const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2063           for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2064             event.push(*itRank,1,msg);
2065           contextClientTmp->sendEvent(event);
2066         }
2067         else contextClientTmp->sendEvent(event);
2068     }
2069   }
2070   CATCH_DUMP_ATTR
2071
2072   //! Server side: Receive a message of client annoucing calendar update
2073   void CContext::recvUpdateCalendar(CEventServer& event)
2074   TRY
2075   {
2076      CBufferIn* buffer=event.subEvents.begin()->buffer;
2077      getCurrent()->recvUpdateCalendar(*buffer);
2078   }
2079   CATCH
2080
2081   //! Server side: Receive a message of client annoucing calendar update
2082   void CContext::recvUpdateCalendar(CBufferIn& buffer)
2083   TRY
2084   {
2085      int step;
2086      buffer>>step;
2087      updateCalendar(step);
2088      if (serviceType_==CServicesManager::GATHERER)
2089      {       
2090        sendUpdateCalendar(step);
2091      }
2092   }
2093   CATCH_DUMP_ATTR
2094
2095   //! Client side: Send a message to create header part of netcdf file
2096   void CContext::sendCreateFileHeader(void)
2097   TRY
2098   {
2099     int nbSrvPools ;
2100     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2101     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2102     else nbSrvPools = 0 ;
2103     CContextClient* contextClientTmp ;
2104
2105     for (int i = 0; i < nbSrvPools; ++i)
2106     {
2107       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2108       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2109       CEventClient event(getType(),EVENT_ID_CREATE_FILE_HEADER);
2110
2111       if (contextClientTmp->isServerLeader())
2112       {
2113         CMessage msg;
2114         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2115         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2116           event.push(*itRank,1,msg) ;
2117         contextClientTmp->sendEvent(event);
2118       }
2119       else contextClientTmp->sendEvent(event);
2120     }
2121   }
2122   CATCH_DUMP_ATTR
2123
2124   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2125   void CContext::recvCreateFileHeader(CEventServer& event)
2126   TRY
2127   {
2128      CBufferIn* buffer=event.subEvents.begin()->buffer;
2129      getCurrent()->recvCreateFileHeader(*buffer);
2130   }
2131   CATCH
2132
2133   //! Server side: Receive a message of client annoucing the creation of header part of netcdf file
2134   void CContext::recvCreateFileHeader(CBufferIn& buffer)
2135   TRY
2136   {
2137      if (serviceType_==CServicesManager::IO_SERVER || serviceType_==CServicesManager::OUT_SERVER) 
2138        createFileHeader();
2139   }
2140   CATCH_DUMP_ATTR
2141
2142   //! Client side: Send a message to do some post processing on server
2143   void CContext::sendProcessingGridOfEnabledFields()
2144   TRY
2145   {
2146     int nbSrvPools ;
2147     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2148     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2149     else nbSrvPools = 0 ;
2150     CContextClient* contextClientTmp ;
2151
2152     for (int i = 0; i < nbSrvPools; ++i)
2153     {
2154       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2155       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2156
2157       CEventClient event(getType(),EVENT_ID_PROCESS_GRID_ENABLED_FIELDS);
2158
2159       if (contextClientTmp->isServerLeader())
2160       {
2161         CMessage msg;
2162         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2163         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2164           event.push(*itRank,1,msg);
2165         contextClientTmp->sendEvent(event);
2166       }
2167       else contextClientTmp->sendEvent(event);
2168     }
2169   }
2170   CATCH_DUMP_ATTR
2171
2172   //! Server side: Receive a message to do some post processing
2173   void CContext::recvProcessingGridOfEnabledFields(CEventServer& event)
2174   TRY
2175   {
2176      CBufferIn* buffer=event.subEvents.begin()->buffer;
2177      // nothing to do, no call ??!!   
2178   }
2179   CATCH
2180
2181   //! Client side: Send a message to do some post processing on server
2182   void CContext::sendPostProcessing()
2183   TRY
2184   {
2185     int nbSrvPools ;
2186     if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2187     else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2188     else nbSrvPools = 0 ;
2189     CContextClient* contextClientTmp ;
2190
2191     for (int i = 0; i < nbSrvPools; ++i)
2192     {
2193       if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2194       else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2195       CEventClient event(getType(),EVENT_ID_POST_PROCESS);
2196       if (contextClientTmp->isServerLeader())
2197       {
2198         CMessage msg;
2199         const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2200         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2201         event.push(*itRank,1,msg);
2202         contextClientTmp->sendEvent(event);
2203       }
2204       else contextClientTmp->sendEvent(event);
2205     }
2206   }
2207   CATCH_DUMP_ATTR
2208
2209   //! Server side: Receive a message to do some post processing
2210   void CContext::recvPostProcessing(CEventServer& event)
2211   TRY
2212   {
2213      CBufferIn* buffer=event.subEvents.begin()->buffer;
2214      getCurrent()->recvPostProcessing(*buffer);
2215   }
2216   CATCH
2217
2218   //! Server side: Receive a message to do some post processing
2219   void CContext::recvPostProcessing(CBufferIn& buffer)
2220   TRY
2221   {
2222      CCalendarWrapper::get(CCalendarWrapper::GetDefName())->createCalendar();
2223      postProcessing();
2224   }
2225   CATCH_DUMP_ATTR
2226
2227 
2228   /*!
2229   \brief Do some simple post processings after parsing xml file
2230      After the xml file (iodef.xml) is parsed, it is necessary to build all relations among
2231   created object, e.g: inhertance among fields, domain, axis. After that, all fiels as well as their parents (reference fields),
2232   which will be written out into netcdf files, are processed
2233   */
2234   void CContext::postProcessing()
2235   TRY
2236   {
2237     if (isPostProcessed) return;
2238
2239      // Make sure the calendar was correctly created
2240      if (!calendar)
2241        ERROR("CContext::postProcessing()", << "A calendar must be defined for the context \"" << getId() << "!\"")
2242      else if (calendar->getTimeStep() == NoneDu)
2243        ERROR("CContext::postProcessing()", << "A timestep must be defined for the context \"" << getId() << "!\"")
2244      // Calendar first update to set the current date equals to the start date
2245      calendar->update(0);
2246
2247      // Find all inheritance in xml structure
2248      this->solveAllInheritance();
2249
2250//      ShowTree(info(10));
2251
2252      // Check if some axis, domains or grids are eligible to for compressed indexed output.
2253      // Warning: This must be done after solving the inheritance and before the rest of post-processing
2254      checkAxisDomainsGridsEligibilityForCompressedOutput();      // only for field written on IO_SERVER service ????
2255
2256      // Check if some automatic time series should be generated
2257      // Warning: This must be done after solving the inheritance and before the rest of post-processing     
2258      prepareTimeseries();
2259
2260      //Initialisation du vecteur 'enabledFiles' contenant la liste des fichiers à sortir.
2261      findEnabledFiles();
2262      findEnabledWriteModeFiles();
2263      findEnabledReadModeFiles();
2264      findEnabledCouplerIn();
2265      findEnabledCouplerOut();
2266      createCouplerInterCommunicator() ;
2267
2268      // Find all enabled fields of each file     
2269      const vector<CField*>&& fileOutField = findAllEnabledFieldsInFiles(this->enabledWriteModeFiles);
2270      const vector<CField*>&& fileInField = findAllEnabledFieldsInFiles(this->enabledReadModeFiles);
2271      const vector<CField*>&& CouplerOutField = findAllEnabledFieldsCouplerOut(this->enabledCouplerOut);
2272      const vector<CField*>&& CouplerInField = findAllEnabledFieldsCouplerIn(this->enabledCouplerIn);
2273
2274
2275
2276      // For now, only read files with client and only one level server
2277      // if (hasClient && !hasServer) findEnabledReadModeFiles();     
2278
2279
2280      // For now, only read files with client and only one level server
2281      // if (hasClient && !hasServer)
2282      //   findAllEnabledFieldsInFiles(this->enabledReadModeFiles);     
2283
2284      if (serviceType_==CServicesManager::CLIENT)
2285      {
2286        initReadFiles();
2287        // Try to read attributes of fields in file then fill in corresponding grid (or domain, axis)
2288        this->readAttributesOfEnabledFieldsInReadModeFiles();
2289      }
2290
2291      // Only search and rebuild all reference objects of enable fields, don't transform
2292      this->solveOnlyRefOfEnabledFields();
2293
2294      // Search and rebuild all reference object of enabled fields, and transform
2295      this->solveAllRefOfEnabledFieldsAndTransform();
2296
2297      // Find all fields with read access from the public API
2298      if (serviceType_==CServicesManager::CLIENT) findFieldsWithReadAccess();
2299      // and solve the all reference for them
2300      if (serviceType_==CServicesManager::CLIENT) solveAllRefOfFieldsWithReadAccess();
2301
2302      isPostProcessed = true;
2303   }
2304   CATCH_DUMP_ATTR
2305
2306   void CContext::createCouplerInterCommunicator(void)
2307   TRY
2308   {
2309      // juste for test now, in future need an scheduler to avoid dead-lock
2310      for(auto it=enabledCouplerOut.begin();it!=enabledCouplerOut.end();++it)
2311      {
2312        (*it)->createInterCommunicator() ;
2313      }
2314
2315      for(auto it=enabledCouplerIn.begin();it!=enabledCouplerIn.end();++it)
2316      {
2317        (*it)->createInterCommunicator() ;
2318      }
2319   }
2320   CATCH_DUMP_ATTR
2321
2322 
2323     //! Client side: Send infomation of active files (files are enabled to write out)
2324   void CContext::sendEnabledFiles(const std::vector<CFile*>& activeFiles)
2325   TRY
2326   {
2327     int size = activeFiles.size();
2328
2329     // In a context, each type has a root definition, e.g: axis, domain, field.
2330     // Every object must be a child of one of these root definition. In this case
2331     // all new file objects created on server must be children of the root "file_definition"
2332     StdString fileDefRoot("file_definition");
2333     CFileGroup* cfgrpPtr = CFileGroup::get(fileDefRoot);
2334
2335     for (int i = 0; i < size; ++i)
2336     {
2337       CFile* f = activeFiles[i];
2338       cfgrpPtr->sendCreateChild(f->getId(),f->getContextClient());
2339       f->sendAllAttributesToServer(f->getContextClient());
2340       f->sendAddAllVariables(f->getContextClient());
2341     }
2342   }
2343   CATCH_DUMP_ATTR
2344
2345   //! Client side: Send information of active fields (ones are written onto files)
2346   void CContext::sendEnabledFieldsInFiles(const std::vector<CFile*>& activeFiles)
2347   TRY
2348   {
2349     int size = activeFiles.size();
2350     for (int i = 0; i < size; ++i)
2351     {
2352       activeFiles[i]->sendEnabledFields(activeFiles[i]->getContextClient());
2353     }
2354   }
2355   CATCH_DUMP_ATTR
2356
2357   //! Client side: Check if the defined axis, domains and grids are eligible for compressed indexed output
2358   void CContext::checkAxisDomainsGridsEligibilityForCompressedOutput()
2359   TRY
2360   {
2361     if (!(serviceType_==CServicesManager::CLIENT || serviceType_==CServicesManager::GATHERER)) return;
2362
2363     const vector<CAxis*> allAxis = CAxis::getAll();
2364     for (vector<CAxis*>::const_iterator it = allAxis.begin(); it != allAxis.end(); it++)
2365       (*it)->checkEligibilityForCompressedOutput();
2366
2367     const vector<CDomain*> allDomains = CDomain::getAll();
2368     for (vector<CDomain*>::const_iterator it = allDomains.begin(); it != allDomains.end(); it++)
2369       (*it)->checkEligibilityForCompressedOutput();
2370
2371     const vector<CGrid*> allGrids = CGrid::getAll();
2372     for (vector<CGrid*>::const_iterator it = allGrids.begin(); it != allGrids.end(); it++)
2373       (*it)->checkEligibilityForCompressedOutput();
2374   }
2375   CATCH_DUMP_ATTR
2376
2377   //! Client side: Prepare the timeseries by adding the necessary files
2378   void CContext::prepareTimeseries()
2379   TRY
2380   {
2381     const std::vector<CFile*> allFiles = CFile::getAll();
2382     for (size_t i = 0; i < allFiles.size(); i++)
2383     {
2384       CFile* file = allFiles[i];
2385
2386       std::vector<CVariable*> fileVars, fieldVars, vars = file->getAllVariables();
2387       for (size_t k = 0; k < vars.size(); k++)
2388       {
2389         CVariable* var = vars[k];
2390
2391         if (var->ts_target.isEmpty()
2392              || var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both)
2393           fileVars.push_back(var);
2394
2395         if (!var->ts_target.isEmpty()
2396              && (var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both))
2397           fieldVars.push_back(var);
2398       }
2399
2400       if (!file->timeseries.isEmpty() && file->timeseries != CFile::timeseries_attr::none)
2401       {
2402         StdString fileNameStr("%file_name%") ;
2403         StdString tsPrefix = !file->ts_prefix.isEmpty() ? file->ts_prefix : fileNameStr ;
2404         
2405         StdString fileName=file->getFileOutputName();
2406         size_t pos=tsPrefix.find(fileNameStr) ;
2407         while (pos!=std::string::npos)
2408         {
2409           tsPrefix=tsPrefix.replace(pos,fileNameStr.size(),fileName) ;
2410           pos=tsPrefix.find(fileNameStr) ;
2411         }
2412       
2413         const std::vector<CField*> allFields = file->getAllFields();
2414         for (size_t j = 0; j < allFields.size(); j++)
2415         {
2416           CField* field = allFields[j];
2417
2418           if (!field->ts_enabled.isEmpty() && field->ts_enabled)
2419           {
2420             CFile* tsFile = CFile::create();
2421             tsFile->duplicateAttributes(file);
2422
2423             // Add variables originating from file and targeted to timeserie file
2424             for (size_t k = 0; k < fileVars.size(); k++)
2425               tsFile->getVirtualVariableGroup()->addChild(fileVars[k]);
2426
2427           
2428             tsFile->name = tsPrefix + "_";
2429             if (!field->name.isEmpty())
2430               tsFile->name.get() += field->name;
2431             else if (field->hasDirectFieldReference()) // We cannot use getBaseFieldReference() just yet
2432               tsFile->name.get() += field->field_ref;
2433             else
2434               tsFile->name.get() += field->getId();
2435
2436             if (!field->ts_split_freq.isEmpty())
2437               tsFile->split_freq = field->ts_split_freq;
2438
2439             CField* tsField = tsFile->addField();
2440             tsField->field_ref = field->getId();
2441
2442             // Add variables originating from file and targeted to timeserie field
2443             for (size_t k = 0; k < fieldVars.size(); k++)
2444               tsField->getVirtualVariableGroup()->addChild(fieldVars[k]);
2445
2446             vars = field->getAllVariables();
2447             for (size_t k = 0; k < vars.size(); k++)
2448             {
2449               CVariable* var = vars[k];
2450
2451               // Add variables originating from field and targeted to timeserie field
2452               if (var->ts_target.isEmpty()
2453                    || var->ts_target == CVariable::ts_target_attr::field || var->ts_target == CVariable::ts_target_attr::both)
2454                 tsField->getVirtualVariableGroup()->addChild(var);
2455
2456               // Add variables originating from field and targeted to timeserie file
2457               if (!var->ts_target.isEmpty()
2458                    && (var->ts_target == CVariable::ts_target_attr::file || var->ts_target == CVariable::ts_target_attr::both))
2459                 tsFile->getVirtualVariableGroup()->addChild(var);
2460             }
2461
2462             tsFile->solveFieldRefInheritance(true);
2463
2464             if (file->timeseries == CFile::timeseries_attr::exclusive)
2465               field->enabled = false;
2466           }
2467         }
2468
2469         // Finally disable the original file is need be
2470         if (file->timeseries == CFile::timeseries_attr::only)
2471          file->enabled = false;
2472       }
2473     }
2474   }
2475   CATCH_DUMP_ATTR
2476
2477   //! Client side: Send information of reference grid of active fields
2478   void CContext::sendRefGrid(const std::vector<CFile*>& activeFiles)
2479   TRY
2480   {
2481     std::set<pair<StdString,CContextClient*>> gridIds;
2482
2483     int sizeFile = activeFiles.size();
2484     CFile* filePtr(NULL);
2485
2486     // Firstly, find all reference grids of all active fields
2487     for (int i = 0; i < sizeFile; ++i)
2488     {
2489       filePtr = activeFiles[i];
2490       std::vector<CField*> enabledFields = filePtr->getEnabledFields();
2491       int sizeField = enabledFields.size();
2492       for (int numField = 0; numField < sizeField; ++numField)
2493       {
2494         if (0 != enabledFields[numField]->getRelGrid())
2495           gridIds.insert(make_pair(CGrid::get(enabledFields[numField]->getRelGrid())->getId(),enabledFields[numField]->getContextClient()));
2496       }
2497     }
2498
2499     // Create all reference grids on server side
2500     StdString gridDefRoot("grid_definition");
2501     CGridGroup* gridPtr = CGridGroup::get(gridDefRoot);
2502     for (auto it = gridIds.begin(); it != gridIds.end(); ++it)
2503     {
2504       gridPtr->sendCreateChild(it->first,it->second);
2505       CGrid::get(it->first)->sendAllAttributesToServer(it->second);
2506       CGrid::get(it->first)->sendAllDomains(it->second);
2507       CGrid::get(it->first)->sendAllAxis(it->second);
2508       CGrid::get(it->first)->sendAllScalars(it->second);
2509     }
2510   }
2511   CATCH_DUMP_ATTR
2512
2513   //! Client side: Send information of reference domain, axis and scalar of active fields
2514   void CContext::sendRefDomainsAxisScalars(const std::vector<CFile*>& activeFiles)
2515   TRY
2516   {
2517     std::set<pair<StdString,CContextClient*>> domainIds, axisIds, scalarIds;
2518
2519     // Find all reference domain and axis of all active fields
2520     int numEnabledFiles = activeFiles.size();
2521     for (int i = 0; i < numEnabledFiles; ++i)
2522     {
2523       std::vector<CField*> enabledFields = activeFiles[i]->getEnabledFields();
2524       int numEnabledFields = enabledFields.size();
2525       for (int j = 0; j < numEnabledFields; ++j)
2526       {
2527         CContextClient* contextClient=enabledFields[j]->getContextClient() ;
2528         const std::vector<StdString>& prDomAxisScalarId = enabledFields[j]->getRefDomainAxisIds();
2529         if ("" != prDomAxisScalarId[0]) domainIds.insert(make_pair(prDomAxisScalarId[0],contextClient));
2530         if ("" != prDomAxisScalarId[1]) axisIds.insert(make_pair(prDomAxisScalarId[1],contextClient));
2531         if ("" != prDomAxisScalarId[2]) scalarIds.insert(make_pair(prDomAxisScalarId[2],contextClient));
2532       }
2533     }
2534
2535     // Create all reference axis on server side
2536     std::set<StdString>::iterator itDom, itAxis, itScalar;
2537     std::set<StdString>::const_iterator itE;
2538
2539     StdString scalarDefRoot("scalar_definition");
2540     CScalarGroup* scalarPtr = CScalarGroup::get(scalarDefRoot);
2541     
2542     for (auto itScalar = scalarIds.begin(); itScalar != scalarIds.end(); ++itScalar)
2543     {
2544       if (!itScalar->first.empty())
2545       {
2546         scalarPtr->sendCreateChild(itScalar->first,itScalar->second);
2547         CScalar::get(itScalar->first)->sendAllAttributesToServer(itScalar->second);
2548       }
2549     }
2550
2551     StdString axiDefRoot("axis_definition");
2552     CAxisGroup* axisPtr = CAxisGroup::get(axiDefRoot);
2553     
2554     for (auto itAxis = axisIds.begin(); itAxis != axisIds.end(); ++itAxis)
2555     {
2556       if (!itAxis->first.empty())
2557       {
2558         axisPtr->sendCreateChild(itAxis->first, itAxis->second);
2559         CAxis::get(itAxis->first)->sendAllAttributesToServer(itAxis->second);
2560       }
2561     }
2562
2563     // Create all reference domains on server side
2564     StdString domDefRoot("domain_definition");
2565     CDomainGroup* domPtr = CDomainGroup::get(domDefRoot);
2566     
2567     for (auto itDom = domainIds.begin(); itDom != domainIds.end(); ++itDom)
2568     {
2569       if (!itDom->first.empty()) {
2570          domPtr->sendCreateChild(itDom->first, itDom->second);
2571          CDomain::get(itDom->first)->sendAllAttributesToServer(itDom->second);
2572       }
2573     }
2574   }
2575   CATCH_DUMP_ATTR
2576
2577   //! Update calendar in each time step
2578   void CContext::updateCalendar(int step)
2579   TRY
2580   {
2581      int prevStep = calendar->getStep();
2582
2583      if (prevStep < step)
2584      {
2585        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2586        {
2587          doPreTimestepOperationsForEnabledReadModeFiles();
2588        }
2589
2590        info(50) << "updateCalendar : before : " << calendar->getCurrentDate() << endl;
2591        calendar->update(step);
2592        info(50) << "updateCalendar : after : " << calendar->getCurrentDate() << endl;
2593  #ifdef XIOS_MEMTRACK_LIGHT
2594        info(50) << " Current memory used by XIOS : "<<  MemTrack::getCurrentMemorySize()*1.0/(1024*1024)<<" Mbyte, at timestep "<<step<<" of context "<<this->getId()<<endl ;
2595  #endif
2596
2597        if (serviceType_==CServicesManager::CLIENT) // For now we only use server level 1 to read data
2598        {
2599          doPostTimestepOperationsForEnabledReadModeFiles();
2600          garbageCollector.invalidate(calendar->getCurrentDate());
2601        }
2602      }
2603      else if (prevStep == step)
2604        info(50) << "updateCalendar: already at step " << step << ", no operation done." << endl;
2605      else // if (prevStep > step)
2606        ERROR("void CContext::updateCalendar(int step)",
2607              << "Illegal calendar update: previous step was " << prevStep << ", new step " << step << "is in the past!")
2608   }
2609   CATCH_DUMP_ATTR
2610
2611   void CContext::initReadFiles(void)
2612   TRY
2613   {
2614      vector<CFile*>::const_iterator it;
2615
2616      for (it=enabledReadModeFiles.begin(); it != enabledReadModeFiles.end(); it++)
2617      {
2618         (*it)->initRead();
2619      }
2620   }
2621   CATCH_DUMP_ATTR
2622
2623   //! Server side: Create header of netcdf file
2624   void CContext::createFileHeader(void)
2625   TRY
2626   {
2627      vector<CFile*>::const_iterator it;
2628
2629      for (it=enabledFiles.begin(); it != enabledFiles.end(); it++)
2630      // for (it=enabledWriteModeFiles.begin(); it != enabledWriteModeFiles.end(); it++)
2631      {
2632         (*it)->initWrite();
2633      }
2634   }
2635   CATCH_DUMP_ATTR
2636
2637   //! Get current context
2638   CContext* CContext::getCurrent(void)
2639   TRY
2640   {
2641     return CObjectFactory::GetObject<CContext>(CObjectFactory::GetCurrentContextId()).get();
2642   }
2643   CATCH
2644
2645   /*!
2646   \brief Set context with an id be the current context
2647   \param [in] id identity of context to be set to current
2648   */
2649   void CContext::setCurrent(const string& id)
2650   TRY
2651   {
2652     CObjectFactory::SetCurrentContextId(id);
2653     CGroupFactory::SetCurrentContextId(id);
2654   }
2655   CATCH
2656
2657  /*!
2658  \brief Create a context with specific id
2659  \param [in] id identity of new context
2660  \return pointer to the new context or already-existed one with identity id
2661  */
2662  CContext* CContext::create(const StdString& id)
2663  TRY
2664  {
2665    CContext::setCurrent(id);
2666
2667    bool hasctxt = CContext::has(id);
2668    CContext* context = CObjectFactory::CreateObject<CContext>(id).get();
2669    getRoot();
2670    if (!hasctxt) CGroupFactory::AddChild(root, context->getShared());
2671
2672#define DECLARE_NODE(Name_, name_) \
2673    C##Name_##Definition::create(C##Name_##Definition::GetDefName());
2674#define DECLARE_NODE_PAR(Name_, name_)
2675#include "node_type.conf"
2676
2677    return (context);
2678  }
2679  CATCH
2680
2681     //! Server side: Receive a message to do some post processing
2682  void CContext::recvRegistry(CEventServer& event)
2683  TRY
2684  {
2685    CBufferIn* buffer=event.subEvents.begin()->buffer;
2686    getCurrent()->recvRegistry(*buffer);
2687  }
2688  CATCH
2689
2690  void CContext::recvRegistry(CBufferIn& buffer)
2691  TRY
2692  {
2693    if (server->intraCommRank==0)
2694    {
2695      CRegistry registry(server->intraComm) ;
2696      registry.fromBuffer(buffer) ;
2697      registryOut->mergeRegistry(registry) ;
2698    }
2699  }
2700  CATCH_DUMP_ATTR
2701
2702  void CContext::sendRegistry(void)
2703  TRY
2704  {
2705    registryOut->hierarchicalGatherRegistry() ;
2706
2707    int nbSrvPools ;
2708    if (serviceType_==CServicesManager::CLIENT) nbSrvPools = 1 ;
2709    else if (serviceType_==CServicesManager::GATHERER) nbSrvPools = this->clientPrimServer.size() ;
2710    else nbSrvPools = 0 ;
2711    CContextClient* contextClientTmp ;
2712
2713    for (int i = 0; i < nbSrvPools; ++i)
2714    {
2715      if (serviceType_==CServicesManager::CLIENT) contextClientTmp = client ;
2716      else if (serviceType_==CServicesManager::GATHERER ) contextClientTmp = clientPrimServer[i] ;
2717
2718      CEventClient event(CContext::GetType(), CContext::EVENT_ID_SEND_REGISTRY);
2719      if (contextClientTmp->isServerLeader())
2720      {
2721        CMessage msg ;
2722        if (contextClientTmp->clientRank==0) msg<<*registryOut ;
2723        const std::list<int>& ranks = contextClientTmp->getRanksServerLeader();
2724        for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2725             event.push(*itRank,1,msg);
2726        contextClientTmp->sendEvent(event);
2727      }
2728      else contextClientTmp->sendEvent(event);
2729    }
2730  }
2731  CATCH_DUMP_ATTR
2732
2733 
2734  void CContext::sendFinalizeClient(CContextClient* contextClient, const string& contextClientId)
2735  TRY
2736  {
2737    CEventClient event(getType(),EVENT_ID_CONTEXT_FINALIZE_CLIENT);
2738    if (contextClient->isServerLeader())
2739    {
2740      CMessage msg;
2741      msg<<contextClientId ;
2742      const std::list<int>& ranks = contextClient->getRanksServerLeader();
2743      for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
2744           event.push(*itRank,1,msg);
2745      contextClient->sendEvent(event);
2746    }
2747    else contextClient->sendEvent(event);
2748  }
2749  CATCH_DUMP_ATTR
2750
2751 
2752  void CContext::recvFinalizeClient(CEventServer& event)
2753  TRY
2754  {
2755    CBufferIn* buffer=event.subEvents.begin()->buffer;
2756    string id;
2757    *buffer>>id;
2758    get(id)->recvFinalizeClient(*buffer);
2759  }
2760  CATCH
2761
2762  void CContext::recvFinalizeClient(CBufferIn& buffer)
2763  TRY
2764  {
2765    countChildContextFinalized_++ ;
2766  }
2767  CATCH_DUMP_ATTR
2768
2769
2770
2771
2772  /*!
2773  * \fn bool CContext::isFinalized(void)
2774  * Context is finalized if it received context post finalize event.
2775  */
2776  bool CContext::isFinalized(void)
2777  TRY
2778  {
2779    return finalized;
2780  }
2781  CATCH_DUMP_ATTR
2782  ///--------------------------------------------------------------
2783  StdString CContext::dumpClassAttributes(void)
2784  {
2785    StdString str;
2786    str.append("enabled files=\"");
2787    int size = this->enabledFiles.size();
2788    for (int i = 0; i < size; ++i)
2789    {
2790      str.append(enabledFiles[i]->getId());
2791      str.append(" ");
2792    }
2793    str.append("\"");
2794    return str;
2795  }
2796
2797} // namespace xios
Note: See TracBrowser for help on using the repository browser.