source: XIOS/dev/dev_ym/XIOS_SERVICES/src/node/context.cpp @ 1761

Last change on this file since 1761 was 1761, checked in by ymipsl, 5 years ago

implementing first guess for service functionnalities.

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