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

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

Some Update on XIOS services
Seems to work on Irène for :

  • first level of servers
  • fisrt + second level of servers
  • attached mode

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