source: XIOS/dev/branch_openmp/src/node/field.cpp @ 1209

Last change on this file since 1209 was 1209, checked in by yushan, 7 years ago

bug corrected. happened when certain threads send 0 elements in the allgatherv call

  • 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
  • Property svn:executable set to *
File size: 47.7 KB
Line 
1#include "field.hpp"
2
3#include "attribute_template.hpp"
4#include "object_template.hpp"
5#include "group_template.hpp"
6
7#include "node_type.hpp"
8#include "calendar_util.hpp"
9#include "message.hpp"
10#include "xios_spl.hpp"
11#include "type.hpp"
12#include "timer.hpp"
13#include "context_client.hpp"
14#include "context_server.hpp"
15#include <set>
16#include "garbage_collector.hpp"
17#include "source_filter.hpp"
18#include "store_filter.hpp"
19#include "file_writer_filter.hpp"
20#include "pass_through_filter.hpp"
21#include "filter_expr_node.hpp"
22#include "lex_parser.hpp"
23#include "temporal_filter.hpp"
24#include "spatial_transform_filter.hpp"
25
26namespace xios{
27
28   /// ////////////////////// Dfinitions ////////////////////// ///
29
30   CField::CField(void)
31      : CObjectTemplate<CField>(), CFieldAttributes()
32      , grid(), file()
33      , written(false)
34      , nstep(0), nstepMax(0)
35      , hasOutputFile(false)
36      , domAxisScalarIds_(vector<StdString>(3,"")), areAllReferenceSolved(false), isReferenceSolved(false)
37      , useCompressedOutput(false)
38      , hasTimeInstant(false)
39      , hasTimeCentered(false)
40      , wasDataAlreadyReceivedFromServer(false)
41      , isEOF(false)
42   { setVirtualVariableGroup(CVariableGroup::create(getId() + "_virtual_variable_group")); }
43
44   CField::CField(const StdString& id)
45      : CObjectTemplate<CField>(id), CFieldAttributes()
46      , grid(), file()
47      , written(false)
48      , nstep(0), nstepMax(0)
49      , hasOutputFile(false)
50      , domAxisScalarIds_(vector<StdString>(3,"")), areAllReferenceSolved(false), isReferenceSolved(false)
51      , useCompressedOutput(false)
52      , hasTimeInstant(false)
53      , hasTimeCentered(false)
54      , wasDataAlreadyReceivedFromServer(false)
55      , isEOF(false)
56   { setVirtualVariableGroup(CVariableGroup::create(getId() + "_virtual_variable_group")); }
57
58   CField::~CField(void)
59   {}
60
61  //----------------------------------------------------------------
62
63   void CField::setVirtualVariableGroup(CVariableGroup* newVVariableGroup)
64   {
65      this->vVariableGroup = newVVariableGroup;
66   }
67
68   CVariableGroup* CField::getVirtualVariableGroup(void) const
69   {
70      return this->vVariableGroup;
71   }
72
73   std::vector<CVariable*> CField::getAllVariables(void) const
74   {
75      return this->vVariableGroup->getAllChildren();
76   }
77
78   void CField::solveDescInheritance(bool apply, const CAttributeMap* const parent)
79   {
80      SuperClassAttribute::setAttributes(parent, apply);
81      this->getVirtualVariableGroup()->solveDescInheritance(apply, NULL);
82   }
83
84  //----------------------------------------------------------------
85
86  bool CField::dispatchEvent(CEventServer& event)
87  {
88    if (SuperClass::dispatchEvent(event)) return true;
89    else
90    {
91      switch(event.type)
92      {
93        case EVENT_ID_UPDATE_DATA :
94          recvUpdateData(event);
95          return true;
96          break;
97
98        case EVENT_ID_READ_DATA :
99          recvReadDataRequest(event);
100          return true;
101          break;
102
103        case EVENT_ID_READ_DATA_READY :
104          recvReadDataReady(event);
105          return true;
106          break;
107
108        case EVENT_ID_ADD_VARIABLE :
109          recvAddVariable(event);
110          return true;
111          break;
112
113        case EVENT_ID_ADD_VARIABLE_GROUP :
114          recvAddVariableGroup(event);
115          return true;
116          break;
117
118        default :
119          ERROR("bool CField::dispatchEvent(CEventServer& event)", << "Unknown Event");
120          return false;
121      }
122    }
123  }
124
125  void CField::sendUpdateData(const CArray<double,1>& data)
126  {
127    CTimer::get("Field : send data").resume();
128
129    CContext* context = CContext::getCurrent();
130    CContextClient* client = context->client;
131
132    CEventClient event(getType(), EVENT_ID_UPDATE_DATA);
133
134    map<int, CArray<int,1> >::iterator it;
135    list<CMessage> list_msg;
136    list<CArray<double,1> > list_data;
137
138    if (!grid->doGridHaveDataDistributed())
139    {
140       if (client->isServerLeader())
141       {
142          for (it = grid->storeIndex_toSrv.begin(); it != grid->storeIndex_toSrv.end(); it++)
143          {
144            int rank = it->first;
145            CArray<int,1>& index = it->second;
146
147            list_msg.push_back(CMessage());
148            list_data.push_back(CArray<double,1>(index.numElements()));
149
150            CArray<double,1>& data_tmp = list_data.back();
151            for (int n = 0; n < data_tmp.numElements(); n++) data_tmp(n) = data(index(n));
152
153            list_msg.back() << getId() << data_tmp;
154            event.push(rank, 1, list_msg.back());
155          }
156          client->sendEvent(event);
157       } 
158       else client->sendEvent(event);
159    }
160    else
161    {
162      for (it = grid->storeIndex_toSrv.begin(); it != grid->storeIndex_toSrv.end(); it++)
163      {
164        int rank = it->first;
165        CArray<int,1>& index = it->second;
166
167        list_msg.push_back(CMessage());
168        list_data.push_back(CArray<double,1>(index.numElements()));
169
170        CArray<double,1>& data_tmp = list_data.back();
171        for (int n = 0; n < data_tmp.numElements(); n++) data_tmp(n) = data(index(n));
172
173        list_msg.back() << getId() << data_tmp;
174        event.push(rank, grid->nbSenders[rank], list_msg.back());
175      }
176      client->sendEvent(event);
177    }
178
179    CTimer::get("Field : send data").suspend();
180  }
181
182  void CField::recvUpdateData(CEventServer& event)
183  {
184    vector<int> ranks;
185    vector<CBufferIn*> buffers;
186
187    list<CEventServer::SSubEvent>::iterator it;
188    string fieldId;
189    CTimer::get("Field : recv data").resume();
190    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
191    {
192      int rank = it->rank;
193      CBufferIn* buffer = it->buffer;
194      *buffer >> fieldId;
195      ranks.push_back(rank);
196      buffers.push_back(buffer);
197    }
198    get(fieldId)->recvUpdateData(ranks,buffers);
199    CTimer::get("Field : recv data").suspend();
200  }
201
202  void  CField::recvUpdateData(vector<int>& ranks, vector<CBufferIn*>& buffers)
203  {
204    if (data_srv.empty())
205    {
206      for (map<int, CArray<size_t, 1> >::iterator it = grid->outIndexFromClient.begin(); it != grid->outIndexFromClient.end(); ++it)
207      {
208        int rank = it->first;
209        data_srv.insert(std::make_pair(rank, CArray<double,1>(it->second.numElements())));
210        foperation_srv.insert(pair<int,boost::shared_ptr<func::CFunctor> >(rank,boost::shared_ptr<func::CFunctor>(new func::CInstant(data_srv[rank]))));
211      }
212    }
213
214    CContext* context = CContext::getCurrent();
215    const CDate& currDate = context->getCalendar()->getCurrentDate();
216    const CDate opeDate      = last_operation_srv +freq_op + freq_operation_srv - freq_op;
217    const CDate writeDate    = last_Write_srv     + freq_write_srv;
218
219    if (opeDate <= currDate)
220    {
221      for (int n = 0; n < ranks.size(); n++)
222      {
223        CArray<double,1> data_tmp;
224        *buffers[n] >> data_tmp;
225        (*foperation_srv[ranks[n]])(data_tmp);
226      }
227      last_operation_srv = currDate;
228    }
229
230    if (writeDate < (currDate + freq_operation_srv))
231    {
232      for (int n = 0; n < ranks.size(); n++)
233      {
234        this->foperation_srv[ranks[n]]->final();
235      }
236
237      last_Write_srv = writeDate;
238      writeField();
239      lastlast_Write_srv = last_Write_srv;
240    }
241  }
242
243  void CField::writeField(void)
244  {
245    if (!getRelFile()->allDomainEmpty)
246    {
247      if (grid->doGridHaveDataToWrite() || getRelFile()->type == CFile::type_attr::one_file)
248      {
249        getRelFile()->checkFile();
250        this->incrementNStep();
251        getRelFile()->getDataOutput()->writeFieldData(CField::get(this));
252      }
253    }
254  }
255
256  bool CField::sendReadDataRequest(const CDate& tsDataRequested)
257  {
258    CContext* context = CContext::getCurrent();
259    CContextClient* client = context->client;
260
261    lastDataRequestedFromServer = tsDataRequested;
262
263    if (!isEOF) // No need to send the request if we already know we are at EOF
264    {
265      CEventClient event(getType(), EVENT_ID_READ_DATA);
266      if (client->isServerLeader())
267      {
268        CMessage msg;
269        msg << getId();
270        const std::list<int>& ranks = client->getRanksServerLeader();
271        for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
272          event.push(*itRank, 1, msg);
273        client->sendEvent(event);
274      }
275      else client->sendEvent(event);
276    }
277    else
278      serverSourceFilter->signalEndOfStream(tsDataRequested);
279
280    return !isEOF;
281  }
282
283  /*!
284  Send request new data read from file if need be, that is the current data is out-of-date.
285  \return true if and only if some data was requested
286  */
287  bool CField::sendReadDataRequestIfNeeded(void)
288  {
289    const CDate& currentDate = CContext::getCurrent()->getCalendar()->getCurrentDate();
290
291    bool dataRequested = false;
292
293    while (currentDate >= lastDataRequestedFromServer)
294    {
295      #pragma omp critical (_output)
296      {
297        info(20) << "currentDate : " << currentDate << endl ;
298        info(20) << "lastDataRequestedFromServer : " << lastDataRequestedFromServer << endl ;
299        info(20) << "file->output_freq.getValue() : " << file->output_freq.getValue() << endl ;
300        info(20) << "lastDataRequestedFromServer + file->output_freq.getValue() : " << lastDataRequestedFromServer + file->output_freq << endl ;
301      }
302      dataRequested |= sendReadDataRequest(lastDataRequestedFromServer + file->output_freq);
303    }
304
305    return dataRequested;
306  }
307
308  void CField::recvReadDataRequest(CEventServer& event)
309  {
310    CBufferIn* buffer = event.subEvents.begin()->buffer;
311    StdString fieldId;
312    *buffer >> fieldId;
313    get(fieldId)->recvReadDataRequest();
314  }
315
316  void CField::recvReadDataRequest(void)
317  {
318    CContext* context = CContext::getCurrent();
319    CContextClient* client = context->client;
320
321    CEventClient event(getType(), EVENT_ID_READ_DATA_READY);
322    std::list<CMessage> msgs;
323
324    bool hasData = readField();
325
326    map<int, CArray<double,1> >::iterator it;
327    if (!grid->doGridHaveDataDistributed())
328    {
329       if (client->isServerLeader())
330       {
331          if (!data_srv.empty())
332          {
333            it = data_srv.begin();
334            const std::list<int>& ranks = client->getRanksServerLeader();
335            for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
336            {
337              msgs.push_back(CMessage());
338              CMessage& msg = msgs.back();
339              msg << getId();
340              if (hasData)
341                msg << getNStep() - 1 << it->second;
342              else
343                msg << int(-1);
344              event.push(*itRank, 1, msg);
345            }
346          }
347          client->sendEvent(event);
348       } 
349       else 
350       {
351          // if (!data_srv.empty())
352          // {
353          //   it = data_srv.begin();
354          //   const std::list<int>& ranks = client->getRanksServerNotLeader();
355          //   for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
356          //   {
357          //     msgs.push_back(CMessage());
358          //     CMessage& msg = msgs.back();
359          //     msg << getId();
360          //     if (hasData)
361          //       msg << getNStep() - 1 << it->second;
362          //     else
363          //       msg << int(-1);
364          //     event.push(*itRank, 1, msg);
365          //   }
366          // }
367          client->sendEvent(event);
368       }
369    }
370    else
371    {
372      for (it = data_srv.begin(); it != data_srv.end(); it++)
373      {
374        msgs.push_back(CMessage());
375        CMessage& msg = msgs.back();
376        msg << getId();
377        if (hasData)
378          msg << getNStep() - 1 << it->second;
379        else
380          msg << int(-1);
381        event.push(it->first, grid->nbSenders[it->first], msg);
382      }
383      client->sendEvent(event);
384    }
385  }
386
387  bool CField::readField(void)
388  {
389    if (!getRelFile()->allDomainEmpty)
390    {
391      if (grid->doGridHaveDataToWrite() || getRelFile()->type == CFile::type_attr::one_file)
392      {
393        if (data_srv.empty())
394        {
395          for (map<int, CArray<size_t, 1> >::iterator it = grid->outIndexFromClient.begin(); it != grid->outIndexFromClient.end(); ++it)
396            data_srv.insert(std::make_pair(it->first, CArray<double,1>(it->second.numElements())));
397        }
398
399        getRelFile()->checkFile();
400        if (!nstepMax)
401        {
402          nstepMax = getRelFile()->getDataInput()->getFieldNbRecords(CField::get(this));
403        }
404
405        this->incrementNStep();
406
407        if (getNStep() > nstepMax && (getRelFile()->cyclic.isEmpty() || !getRelFile()->cyclic) )
408          return false;
409
410        getRelFile()->getDataInput()->readFieldData(CField::get(this));
411      }
412    }
413
414    return true;
415  }
416
417  void CField::recvReadDataReady(CEventServer& event)
418  {
419    string fieldId;
420    vector<int> ranks;
421    vector<CBufferIn*> buffers;
422
423    list<CEventServer::SSubEvent>::iterator it;
424    for (it = event.subEvents.begin(); it != event.subEvents.end(); ++it)
425    {
426      ranks.push_back(it->rank);
427      CBufferIn* buffer = it->buffer;
428      *buffer >> fieldId;
429      buffers.push_back(buffer);
430    }
431    get(fieldId)->recvReadDataReady(ranks, buffers);
432  }
433
434  void CField::recvReadDataReady(vector<int> ranks, vector<CBufferIn*> buffers)
435  {
436    CContext* context = CContext::getCurrent();
437    int record;
438    std::map<int, CArray<double,1> > data;
439
440    for (int i = 0; i < ranks.size(); i++)
441    {
442      int rank = ranks[i];
443      *buffers[i] >> record;
444      isEOF = (record == int(-1));
445
446      if (!isEOF)
447        *buffers[i] >> data[rank];
448      else
449        break;
450    }
451
452    if (wasDataAlreadyReceivedFromServer)
453      lastDataReceivedFromServer = lastDataReceivedFromServer + file->output_freq;
454    else
455    {
456      lastDataReceivedFromServer = context->getCalendar()->getInitDate();
457      wasDataAlreadyReceivedFromServer = true;
458    }
459
460    if (isEOF)
461      serverSourceFilter->signalEndOfStream(lastDataReceivedFromServer);
462    else
463      serverSourceFilter->streamDataFromServer(lastDataReceivedFromServer, data);
464  }
465
466   //----------------------------------------------------------------
467
468   void CField::setRelFile(CFile* _file)
469   {
470      this->file = _file;
471      hasOutputFile = true;
472   }
473
474   //----------------------------------------------------------------
475
476   StdString CField::GetName(void)    { return StdString("field"); }
477   StdString CField::GetDefName(void) { return CField::GetName(); }
478   ENodeType CField::GetType(void)    { return eField; }
479
480   //----------------------------------------------------------------
481
482   CGrid* CField::getRelGrid(void) const
483   {
484      return this->grid;
485   }
486
487   //----------------------------------------------------------------
488
489   CFile* CField::getRelFile(void) const
490   {
491      return this->file;
492   }
493
494   int CField::getNStep(void) const
495   {
496      return this->nstep;
497   }
498
499   func::CFunctor::ETimeType CField::getOperationTimeType() const
500   {
501     return operationTimeType;
502   }
503
504   //----------------------------------------------------------------
505
506   void CField::incrementNStep(void)
507   {
508      this->nstep++;
509   }
510
511   void CField::resetNStep(int nstep /*= 0*/)
512   {
513      this->nstep = nstep;
514   }
515
516   void CField::resetNStepMax(void)
517   {
518      this->nstepMax = 0;
519   }
520
521   //----------------------------------------------------------------
522
523   bool CField::isActive(bool atCurrentTimestep /*= false*/) const
524   {
525      if (clientSourceFilter)
526        return atCurrentTimestep ? clientSourceFilter->isDataExpected(CContext::getCurrent()->getCalendar()->getCurrentDate()) : true;
527      else if (storeFilter)
528        return true;
529      else if (instantDataFilter)
530        ERROR("bool CField::isActive(bool atCurrentTimestep)",
531              << "Impossible to check if field [ id = " << getId() << " ] is active as it cannot be used to receive nor send data.");
532
533      return false;
534   }
535
536   //----------------------------------------------------------------
537
538   bool CField::wasWritten() const
539   {
540     return written;
541   }
542
543   void CField::setWritten()
544   {
545     written = true;
546   }
547
548   //----------------------------------------------------------------
549
550   bool CField::getUseCompressedOutput() const
551   {
552     return useCompressedOutput;
553   }
554
555   void CField::setUseCompressedOutput()
556   {
557     useCompressedOutput = true;
558   }
559
560   //----------------------------------------------------------------
561
562   boost::shared_ptr<COutputPin> CField::getInstantDataFilter()
563   {
564     return instantDataFilter;
565   }
566
567   //----------------------------------------------------------------
568
569   void CField::solveOnlyReferenceEnabledField(bool doSending2Server)
570   {
571     CContext* context = CContext::getCurrent();
572     if (!isReferenceSolved)
573     {
574        isReferenceSolved = true;
575
576        if (context->hasClient)
577        {
578          solveRefInheritance(true);
579          if (hasDirectFieldReference()) getDirectFieldReference()->solveOnlyReferenceEnabledField(false);
580        }
581        else if (context->hasServer)
582          solveServerOperation();
583
584        solveGridReference();
585
586       if (context->hasClient)
587       {
588         solveGenerateGrid();
589         buildGridTransformationGraph();
590       }
591     }
592   }
593
594   /*!
595     Build up graph of grids which plays role of destination and source in grid transformation
596     This function should be called before \func solveGridReference()
597   */
598   void CField::buildGridTransformationGraph()
599   {
600     CContext* context = CContext::getCurrent();
601     if (context->hasClient)
602     {
603       if (grid && !grid->isTransformed() && hasDirectFieldReference() && grid != getDirectFieldReference()->grid)
604       {
605         grid->addTransGridSource(getDirectFieldReference()->grid);
606       }
607     }
608   }
609
610   /*!
611     Generate a new grid destination if there are more than one grid source pointing to a same grid destination
612   */
613   void CField::generateNewTransformationGridDest()
614   {
615     CContext* context = CContext::getCurrent();
616     if (context->hasClient)
617     {
618       std::map<CGrid*,std::pair<bool,StdString> >& gridSrcMap = grid->getTransGridSource();
619       if (1 < gridSrcMap.size())
620       {
621         // Search for grid source
622         CGrid* gridSrc = grid;
623         CField* currField = this;
624         std::vector<CField*> hieraField;
625
626         while (currField->hasDirectFieldReference() && (gridSrc == grid))
627         {
628           hieraField.push_back(currField);
629           CField* tmp = currField->getDirectFieldReference();
630           currField = tmp;
631           gridSrc = currField->grid;
632         }
633
634         if (gridSrcMap.end() != gridSrcMap.find(gridSrc))
635         {
636           CGrid* gridTmp;
637           std::pair<bool,StdString> newGridDest = gridSrcMap[gridSrc];
638           if (newGridDest.first)
639           {
640             StdString newIdGridDest = newGridDest.second;
641             if (!CGrid::has(newIdGridDest))
642             {
643                ERROR("CGrid* CGrid::generateNewTransformationGridDest()",
644                  << " Something wrong happened! Grid whose id " << newIdGridDest
645                  << "should exist ");
646             }
647             gridTmp = CGrid::get(newIdGridDest);
648           }
649           else
650           {
651             StdString newIdGridDest = CGrid::generateId(gridSrc, grid);
652             gridTmp = CGrid::cloneGrid(newIdGridDest, grid);
653
654             (gridSrcMap[gridSrc]).first = true;
655             (gridSrcMap[gridSrc]).second = newIdGridDest;
656           }
657
658           // Update all descendants
659           for (std::vector<CField*>::iterator it = hieraField.begin(); it != hieraField.end(); ++it)
660           {
661             (*it)->grid = gridTmp;
662             (*it)->updateRef((*it)->grid);
663           }
664         }
665       }
666     }
667   }
668
669   void CField::updateRef(CGrid* grid)
670   {
671     if (!grid_ref.isEmpty()) grid_ref.setValue(grid->getId());
672     else
673     {
674       std::vector<CAxis*> axisTmp = grid->getAxis();
675       std::vector<CDomain*> domainTmp = grid->getDomains();
676       if ((1<axisTmp.size()) || (1<domainTmp.size()))
677         ERROR("void CField::updateRef(CGrid* grid)",
678           << "More than one domain or axis is available for domain_ref/axis_ref of field " << this->getId());
679
680       if ((!domain_ref.isEmpty()) && (domainTmp.empty()))
681         ERROR("void CField::updateRef(CGrid* grid)",
682           << "Incoherent between available domain and domain_ref of field " << this->getId());
683       if ((!axis_ref.isEmpty()) && (axisTmp.empty()))
684         ERROR("void CField::updateRef(CGrid* grid)",
685           << "Incoherent between available axis and axis_ref of field " << this->getId());
686
687       if (!domain_ref.isEmpty()) domain_ref.setValue(domainTmp[0]->getId());
688       if (!axis_ref.isEmpty()) axis_ref.setValue(axisTmp[0]->getId());
689     }
690   }
691
692   void CField::solveAllReferenceEnabledField(bool doSending2Server)
693   {
694     CContext* context = CContext::getCurrent();
695     solveOnlyReferenceEnabledField(doSending2Server);
696     int myRank;
697     MPI_Comm_rank(context->client->intraComm, &myRank);
698
699     if (!areAllReferenceSolved)
700     {
701        areAllReferenceSolved = true;
702
703        if (context->hasClient)
704        {
705          solveRefInheritance(true);
706          if (hasDirectFieldReference()) getDirectFieldReference()->solveAllReferenceEnabledField(false);
707        }
708        else if (context->hasServer)
709          solveServerOperation();
710
711        solveGridReference();
712     }
713
714     solveGridDomainAxisRef(doSending2Server);
715
716     if (context->hasClient)
717     {
718       solveTransformedGrid();
719     }
720
721     solveCheckMaskIndex(doSending2Server);
722   }
723
724   std::map<int, StdSize> CField::getGridAttributesBufferSize()
725   {
726     return grid->getAttributesBufferSize();
727   }
728
729   std::map<int, StdSize> CField::getGridDataBufferSize()
730   {
731     return grid->getDataBufferSize(getId());
732   }
733
734   //----------------------------------------------------------------
735
736   void CField::solveServerOperation(void)
737   {
738      CContext* context = CContext::getCurrent();
739
740      if (!context->hasServer || !hasOutputFile) return;
741
742      if (freq_op.isEmpty())
743        freq_op.setValue(TimeStep);
744
745      if (freq_offset.isEmpty())
746        freq_offset.setValue(NoneDu);
747
748      freq_operation_srv = file->output_freq.getValue();
749      freq_write_srv     = file->output_freq.getValue();
750
751      lastlast_Write_srv = context->getCalendar()->getInitDate();
752      last_Write_srv     = context->getCalendar()->getInitDate();
753      last_operation_srv = context->getCalendar()->getInitDate();
754
755      const CDuration toffset = freq_operation_srv - freq_offset.getValue() - context->getCalendar()->getTimeStep();
756      last_operation_srv     = last_operation_srv - toffset;
757
758      if (operation.isEmpty())
759        ERROR("void CField::solveServerOperation(void)",
760              << "An operation must be defined for field \"" << getId() << "\".");
761
762      boost::shared_ptr<func::CFunctor> functor;
763      CArray<double, 1> dummyData;
764
765#define DECLARE_FUNCTOR(MType, mtype) \
766      if (operation.getValue().compare(#mtype) == 0) \
767      { \
768        functor.reset(new func::C##MType(dummyData)); \
769      }
770
771#include "functor_type.conf"
772
773      if (!functor)
774        ERROR("void CField::solveServerOperation(void)",
775              << "\"" << operation << "\" is not a valid operation.");
776
777      operationTimeType = functor->timeType();
778   }
779
780   //----------------------------------------------------------------
781
782   /*!
783    * Constructs the graph filter for the field, enabling or not the data output.
784    * This method should not be called more than once with enableOutput equal to true.
785    *
786    * \param gc the garbage collector to use when building the filter graph
787    * \param enableOutput must be true when the field data is to be
788    *                     read by the client or/and written to a file
789    */
790   void CField::buildFilterGraph(CGarbageCollector& gc, bool enableOutput)
791   {
792     if (!areAllReferenceSolved) solveAllReferenceEnabledField(false);
793
794     const bool detectMissingValues = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
795     const double defaultValue  = detectMissingValues ? default_value : (!default_value.isEmpty() ? default_value : 0.0);
796
797     // Start by building a filter which can provide the field's instant data
798     if (!instantDataFilter)
799     {
800       // Check if we have an expression to parse
801       if (hasExpression())
802       {
803         boost::scoped_ptr<IFilterExprNode> expr(parseExpr(getExpression() + '\0'));
804         boost::shared_ptr<COutputPin> filter = expr->reduce(gc, *this);
805
806         // Check if a spatial transformation is needed
807         if (!field_ref.isEmpty())
808         {
809           CGrid* gridRef = CField::get(field_ref)->grid;
810
811           if (grid && grid != gridRef && grid->hasTransform())
812           {
813             std::pair<boost::shared_ptr<CFilter>, boost::shared_ptr<CFilter> > filters = CSpatialTransformFilter::buildFilterGraph(gc, gridRef, grid, detectMissingValues, defaultValue);
814
815             filter->connectOutput(filters.first, 0);
816             filter = filters.second;
817           }
818         }
819
820         instantDataFilter = filter;
821       }
822       // Check if we have a reference on another field
823       else if (!field_ref.isEmpty())
824         instantDataFilter = getFieldReference(gc);
825       // Check if the data is to be read from a file
826       else if (file && !file->mode.isEmpty() && file->mode == CFile::mode_attr::read)
827         instantDataFilter = serverSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid,
828                                                                                                     freq_offset.isEmpty() ? NoneDu : freq_offset,
829                                                                                                     true,
830                                                                                                     detectMissingValues, defaultValue));
831       else // The data might be passed from the model
832       {
833          if (check_if_active.isEmpty()) check_if_active = false;
834          instantDataFilter = clientSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid, NoneDu, false,
835                                                                                                      detectMissingValues, defaultValue));
836       }
837     }
838
839     // If the field data is to be read by the client or/and written to a file
840     if (enableOutput && !storeFilter && !fileWriterFilter)
841     {
842       if (!read_access.isEmpty() && read_access)
843       {
844         storeFilter = boost::shared_ptr<CStoreFilter>(new CStoreFilter(gc, CContext::getCurrent(), grid,
845                                                                        detectMissingValues, defaultValue));
846         instantDataFilter->connectOutput(storeFilter, 0);
847       }
848
849       if (file && (file->mode.isEmpty() || file->mode == CFile::mode_attr::write))
850       {
851         fileWriterFilter = boost::shared_ptr<CFileWriterFilter>(new CFileWriterFilter(gc, this));
852         getTemporalDataFilter(gc, file->output_freq)->connectOutput(fileWriterFilter, 0);
853       }
854     }
855   }
856
857   /*!
858    * Returns the filter needed to handle the field reference.
859    * This method should only be called when building the filter graph corresponding to the field.
860    *
861    * \param gc the garbage collector to use
862    * \return the output pin corresponding to the field reference
863    */
864   boost::shared_ptr<COutputPin> CField::getFieldReference(CGarbageCollector& gc)
865   {
866     if (instantDataFilter || field_ref.isEmpty())
867       ERROR("COutputPin* CField::getFieldReference(CGarbageCollector& gc)",
868             "Impossible to get the field reference for a field which has already been parsed or which does not have a field_ref.");
869
870     CField* fieldRef = CField::get(field_ref);
871     fieldRef->buildFilterGraph(gc, false);
872
873     std::pair<boost::shared_ptr<CFilter>, boost::shared_ptr<CFilter> > filters;
874     // Check if a spatial transformation is needed
875     if (grid && grid != fieldRef->grid && grid->hasTransform())
876     {       
877       bool hasMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
878       double defaultValue  = hasMissingValue ? default_value : (!default_value.isEmpty() ? default_value : 0.0);                               
879       filters = CSpatialTransformFilter::buildFilterGraph(gc, fieldRef->grid, grid, hasMissingValue, defaultValue);
880     }
881     else
882       filters.first = filters.second = boost::shared_ptr<CFilter>(new CPassThroughFilter(gc));
883
884     fieldRef->getInstantDataFilter()->connectOutput(filters.first, 0);
885
886     return filters.second;
887   }
888
889   /*!
890    * Returns the filter needed to handle a self reference in the field's expression.
891    * If the needed filter does not exist, it is created, otherwise it is reused.
892    * This method should only be called when building the filter graph corresponding
893    * to the field's expression.
894    *
895    * \param gc the garbage collector to use
896    * \return the output pin corresponding to a self reference
897    */
898   boost::shared_ptr<COutputPin> CField::getSelfReference(CGarbageCollector& gc)
899   {
900     if (instantDataFilter || !hasExpression())
901       ERROR("COutputPin* CField::getSelfReference(CGarbageCollector& gc)",
902             "Impossible to add a self reference to a field which has already been parsed or which does not have an expression.");
903
904     if (!selfReferenceFilter)
905     {
906       const bool detectMissingValues = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
907       const double defaultValue  = detectMissingValues ? default_value : (!default_value.isEmpty() ? default_value : 0.0);
908
909       if (file && !file->mode.isEmpty() && file->mode == CFile::mode_attr::read)
910       {
911         if (!serverSourceFilter)
912           serverSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid,
913                                                                                   freq_offset.isEmpty() ? NoneDu : freq_offset,
914                                                                                   true,
915                                                                                   detectMissingValues, defaultValue));
916
917         selfReferenceFilter = serverSourceFilter;
918       }
919       else if (!field_ref.isEmpty())
920       {
921         CField* fieldRef = CField::get(field_ref);
922         fieldRef->buildFilterGraph(gc, false); 
923         selfReferenceFilter = fieldRef->getInstantDataFilter();
924       }
925       else
926       {
927         if (!clientSourceFilter)
928         {
929           if (check_if_active.isEmpty()) check_if_active = false;
930           clientSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid, NoneDu, false,
931                                                                                   detectMissingValues, defaultValue));
932         }
933
934         selfReferenceFilter = clientSourceFilter;
935       }
936     }
937
938     return selfReferenceFilter;
939   }
940
941   /*!
942    * Returns the temporal filter corresponding to the field's temporal operation
943    * for the specified operation frequency. The filter is created if it does not
944    * exist, otherwise it is reused.
945    *
946    * \param gc the garbage collector to use
947    * \param outFreq the operation frequency, i.e. the frequency at which the output data will be computed
948    * \return the output pin corresponding to the requested temporal filter
949    */
950   boost::shared_ptr<COutputPin> CField::getTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)
951   {
952     std::map<CDuration, boost::shared_ptr<COutputPin> >::iterator it = temporalDataFilters.find(outFreq);
953
954     if (it == temporalDataFilters.end())
955     {
956       if (operation.isEmpty())
957         ERROR("void CField::getTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)",
958               << "An operation must be defined for field \"" << getId() << "\".");
959
960       if (freq_op.isEmpty())
961         freq_op.setValue(TimeStep);
962       if (freq_offset.isEmpty())
963         freq_offset.setValue(NoneDu);
964
965       const bool detectMissingValues = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
966       
967       boost::shared_ptr<CTemporalFilter> temporalFilter(new CTemporalFilter(gc, operation,
968                                                                             CContext::getCurrent()->getCalendar()->getInitDate(),
969                                                                             freq_op, freq_offset, outFreq,
970                                                                             detectMissingValues, detectMissingValues ? default_value : 0.0));
971       instantDataFilter->connectOutput(temporalFilter, 0);
972
973       it = temporalDataFilters.insert(std::make_pair(outFreq, temporalFilter)).first;
974     }
975
976     return it->second;
977   }
978
979  /*!
980    * Returns the temporal filter corresponding to the field's temporal operation
981    * for the specified operation frequency.
982    *
983    * \param gc the garbage collector to use
984    * \param outFreq the operation frequency, i.e. the frequency at which the output data will be computed
985    * \return the output pin corresponding to the requested temporal filter
986    */
987   
988   boost::shared_ptr<COutputPin> CField::getSelfTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)
989   {
990     if (instantDataFilter || !hasExpression())
991       ERROR("COutputPin* CField::getSelfTemporalDataFilter(CGarbageCollector& gc)",
992             "Impossible to add a self reference to a field which has already been parsed or which does not have an expression.");
993
994     if (!selfReferenceFilter) getSelfReference(gc) ;
995
996     if (serverSourceFilter || clientSourceFilter)
997     {
998       if (operation.isEmpty())
999         ERROR("void CField::getSelfTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)",
1000               << "An operation must be defined for field \"" << getId() << "\".");
1001
1002       if (freq_op.isEmpty()) freq_op.setValue(TimeStep);
1003       if (freq_offset.isEmpty()) freq_offset.setValue(NoneDu);
1004
1005       const bool detectMissingValues = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
1006
1007       boost::shared_ptr<CTemporalFilter> temporalFilter(new CTemporalFilter(gc, operation,
1008                                                                             CContext::getCurrent()->getCalendar()->getInitDate(),
1009                                                                             freq_op, freq_offset, outFreq,
1010                                                                             detectMissingValues, detectMissingValues ? default_value : 0.0));
1011       selfReferenceFilter->connectOutput(temporalFilter, 0);
1012       return temporalFilter ;
1013     }
1014     else if (!field_ref.isEmpty())
1015     {
1016       CField* fieldRef = CField::get(field_ref);
1017       fieldRef->buildFilterGraph(gc, false); 
1018       return fieldRef->getTemporalDataFilter(gc, outFreq) ;
1019     }
1020  }
1021
1022   //----------------------------------------------------------------
1023/*
1024   void CField::fromBinary(StdIStream& is)
1025   {
1026      SuperClass::fromBinary(is);
1027#define CLEAR_ATT(name_)\
1028      SuperClassAttribute::operator[](#name_)->reset()
1029
1030         CLEAR_ATT(domain_ref);
1031         CLEAR_ATT(axis_ref);
1032#undef CLEAR_ATT
1033
1034   }
1035*/
1036   //----------------------------------------------------------------
1037
1038   void CField::solveGridReference(void)
1039   {
1040      if (grid_ref.isEmpty() && domain_ref.isEmpty() && axis_ref.isEmpty() && scalar_ref.isEmpty())
1041      {
1042        ERROR("CField::solveGridReference(void)",
1043              << "A grid must be defined for field '" << getFieldOutputName() << "' .");
1044      }
1045      else if (!grid_ref.isEmpty() && (!domain_ref.isEmpty() || !axis_ref.isEmpty() || !scalar_ref.isEmpty()))
1046      {
1047        ERROR("CField::solveGridReference(void)",
1048              << "Field '" << getFieldOutputName() << "' has both a grid and a domain/axis/scalar." << std::endl
1049              << "Please define either 'grid_ref' or 'domain_ref'/'axis_ref'/'scalar_ref'.");
1050      }
1051
1052      if (grid_ref.isEmpty())
1053      {
1054        std::vector<CDomain*> vecDom;
1055        std::vector<CAxis*> vecAxis;
1056        std::vector<CScalar*> vecScalar;
1057        std::vector<int> axisDomainOrderTmp;
1058       
1059        if (!domain_ref.isEmpty())
1060        {
1061          StdString tmp = domain_ref.getValue();
1062          if (CDomain::has(domain_ref))
1063          {
1064            vecDom.push_back(CDomain::get(domain_ref));
1065            axisDomainOrderTmp.push_back(2);
1066          }
1067          else
1068            ERROR("CField::solveGridReference(void)",
1069                  << "Invalid reference to domain '" << domain_ref.getValue() << "'.");
1070        }
1071
1072        if (!axis_ref.isEmpty())
1073        {
1074          if (CAxis::has(axis_ref))
1075          {
1076            vecAxis.push_back(CAxis::get(axis_ref));
1077            axisDomainOrderTmp.push_back(1);
1078          }
1079          else
1080            ERROR("CField::solveGridReference(void)",
1081                  << "Invalid reference to axis '" << axis_ref.getValue() << "'.");
1082        }
1083
1084        if (!scalar_ref.isEmpty())
1085        {
1086          if (CScalar::has(scalar_ref))
1087          {
1088            vecScalar.push_back(CScalar::get(scalar_ref));
1089            axisDomainOrderTmp.push_back(0);
1090          }
1091          else
1092            ERROR("CField::solveGridReference(void)",
1093                  << "Invalid reference to scalar '" << scalar_ref.getValue() << "'.");
1094        }
1095       
1096        CArray<int,1> axisDomainOrder(axisDomainOrderTmp.size());
1097        for (int idx = 0; idx < axisDomainOrderTmp.size(); ++idx)
1098        {
1099          axisDomainOrder(idx) = axisDomainOrderTmp[idx];
1100        }
1101
1102        // Warning: the gridId shouldn't be set as the grid_ref since it could be inherited
1103        StdString gridId = CGrid::generateId(vecDom, vecAxis, vecScalar,axisDomainOrder);
1104        if (CGrid::has(gridId))
1105          this->grid = CGrid::get(gridId);
1106        else
1107          this->grid = CGrid::createGrid(gridId, vecDom, vecAxis, vecScalar,axisDomainOrder);
1108      }
1109      else
1110      {
1111        if (CGrid::has(grid_ref))
1112          this->grid = CGrid::get(grid_ref);
1113        else
1114          ERROR("CField::solveGridReference(void)",
1115                << "Invalid reference to grid '" << grid_ref.getValue() << "'.");
1116      }
1117   }
1118
1119   void CField::solveGridDomainAxisRef(bool checkAtt)
1120   {
1121     grid->solveDomainAxisRef(checkAtt);
1122   }
1123
1124   void CField::solveCheckMaskIndex(bool doSendingIndex)
1125   {
1126     grid->checkMaskIndex(doSendingIndex);
1127   }
1128
1129   void CField::solveTransformedGrid()
1130   {
1131     if (grid && !grid->isTransformed() && hasDirectFieldReference() && grid != getDirectFieldReference()->grid)
1132     {
1133       std::vector<CGrid*> grids;
1134       // Source grid
1135       grids.push_back(getDirectFieldReference()->grid);
1136       // Intermediate grids
1137       if (!grid_path.isEmpty())
1138       {
1139         std::string gridId;
1140         size_t start = 0, end;
1141
1142         do
1143         {
1144           end = grid_path.getValue().find(',', start);
1145           if (end != std::string::npos)
1146           {
1147             gridId = grid_path.getValue().substr(start, end - start);
1148             start = end + 1;
1149           }
1150           else
1151             gridId = grid_path.getValue().substr(start);
1152
1153           if (!CGrid::has(gridId))
1154             ERROR("void CField::solveTransformedGrid()",
1155                   << "Invalid grid_path, the grid '" << gridId << "' does not exist.");
1156
1157           grids.push_back(CGrid::get(gridId));
1158         }
1159         while (end != std::string::npos);
1160       }
1161       // Destination grid
1162       grids.push_back(grid);
1163
1164       for (size_t i = 0, count = grids.size() - 1; i < count; ++i)
1165       {
1166         CGrid *gridSrc  = grids[i];
1167         CGrid *gridDest = grids[i + 1];
1168         if (!gridDest->isTransformed())
1169           gridDest->transformGrid(gridSrc);
1170       }
1171     }
1172     else if (grid && grid->hasTransform() && !grid->isTransformed())
1173     {
1174       // Temporarily deactivate the self-transformation of grid
1175       //grid->transformGrid(grid);
1176     }
1177   }
1178
1179   void CField::solveGenerateGrid()
1180   {
1181     if (grid && !grid->isTransformed() && hasDirectFieldReference() && grid != getDirectFieldReference()->grid)
1182       grid->completeGrid(getDirectFieldReference()->grid);
1183     else
1184       grid->completeGrid();
1185   }
1186
1187   void CField::solveGridDomainAxisBaseRef()
1188   {
1189     grid->solveDomainAxisRef(false);
1190     grid->solveDomainAxisBaseRef();
1191   }
1192
1193   ///-------------------------------------------------------------------
1194
1195   template <>
1196   void CGroupTemplate<CField, CFieldGroup, CFieldAttributes>::solveRefInheritance(void)
1197   {
1198      if (this->group_ref.isEmpty()) return;
1199      StdString gref = this->group_ref.getValue();
1200
1201      if (!CFieldGroup::has(gref))
1202         ERROR("CGroupTemplate<CField, CFieldGroup, CFieldAttributes>::solveRefInheritance(void)",
1203               << "[ gref = " << gref << "]"
1204               << " invalid group name !");
1205
1206      CFieldGroup* group = CFieldGroup::get(gref);
1207      CFieldGroup* owner = CFieldGroup::get(boost::polymorphic_downcast<CFieldGroup*>(this));
1208
1209      std::vector<CField*> allChildren  = group->getAllChildren();
1210      std::vector<CField*>::iterator it = allChildren.begin(), end = allChildren.end();
1211
1212      for (; it != end; it++)
1213      {
1214         CField* child = *it;
1215         if (child->hasId()) owner->createChild()->field_ref.setValue(child->getId());
1216
1217      }
1218   }
1219
1220   void CField::scaleFactorAddOffset(double scaleFactor, double addOffset)
1221   {
1222     map<int, CArray<double,1> >::iterator it;
1223     for (it = data_srv.begin(); it != data_srv.end(); it++) it->second = (it->second - addOffset) / scaleFactor;
1224   }
1225
1226   void CField::invertScaleFactorAddOffset(double scaleFactor, double addOffset)
1227   {
1228     map<int, CArray<double,1> >::iterator it;
1229     for (it = data_srv.begin(); it != data_srv.end(); it++) it->second = it->second * scaleFactor + addOffset;
1230   }
1231
1232   void CField::outputField(CArray<double,3>& fieldOut)
1233   {
1234      map<int, CArray<double,1> >::iterator it;
1235      for (it = data_srv.begin(); it != data_srv.end(); it++)
1236      {
1237        grid->outputField(it->first, it->second, fieldOut.dataFirst());
1238      }
1239   }
1240
1241   void CField::outputField(CArray<double,2>& fieldOut)
1242   {
1243      map<int, CArray<double,1> >::iterator it;
1244      for(it=data_srv.begin();it!=data_srv.end();it++)
1245      {
1246         grid->outputField(it->first, it->second, fieldOut.dataFirst());
1247      }
1248   }
1249
1250   void CField::outputField(CArray<double,1>& fieldOut)
1251   {
1252      map<int, CArray<double,1> >::iterator it;
1253
1254      for (it = data_srv.begin(); it != data_srv.end(); it++)
1255      {
1256         grid->outputField(it->first, it->second, fieldOut.dataFirst());
1257      }
1258   }
1259
1260   void CField::inputField(CArray<double,3>& fieldOut)
1261   {
1262      map<int, CArray<double,1> >::iterator it;
1263      for (it = data_srv.begin(); it != data_srv.end(); it++)
1264      {
1265        grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1266      }
1267   }
1268
1269   void CField::inputField(CArray<double,2>& fieldOut)
1270   {
1271      map<int, CArray<double,1> >::iterator it;
1272      for(it = data_srv.begin(); it != data_srv.end(); it++)
1273      {
1274         grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1275      }
1276   }
1277
1278   void CField::inputField(CArray<double,1>& fieldOut)
1279   {
1280      map<int, CArray<double,1> >::iterator it;
1281      for (it = data_srv.begin(); it != data_srv.end(); it++)
1282      {
1283         grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1284      }
1285   }
1286
1287   void CField::outputCompressedField(CArray<double,1>& fieldOut)
1288   {
1289      map<int, CArray<double,1> >::iterator it;
1290
1291      for (it = data_srv.begin(); it != data_srv.end(); it++)
1292      {
1293         grid->outputCompressedField(it->first, it->second, fieldOut.dataFirst());
1294      }
1295   }
1296
1297   ///-------------------------------------------------------------------
1298
1299   void CField::parse(xml::CXMLNode& node)
1300   {
1301      SuperClass::parse(node);
1302      if (!node.getContent(this->content))
1303      {
1304        if (node.goToChildElement())
1305        {
1306          do
1307          {
1308            if (node.getElementName() == "variable" || node.getElementName() == "variable_group") this->getVirtualVariableGroup()->parseChild(node);
1309          } while (node.goToNextElement());
1310          node.goToParentElement();
1311        }
1312      }
1313    }
1314
1315   /*!
1316     This function retrieves Id of corresponding domain_ref and axis_ref (if any)
1317   of a field. In some cases, only domain exists but axis doesn't
1318   \return pair of Domain and Axis id
1319   */
1320   const std::vector<StdString>& CField::getRefDomainAxisIds()
1321   {
1322     CGrid* cgPtr = getRelGrid();
1323     if (NULL != cgPtr)
1324     {
1325       std::vector<StdString>::iterator it;
1326       if (!domain_ref.isEmpty())
1327       {
1328         std::vector<StdString> domainList = cgPtr->getDomainList();
1329         it = std::find(domainList.begin(), domainList.end(), domain_ref.getValue());
1330         if (domainList.end() != it) domAxisScalarIds_[0] = *it;
1331       }
1332
1333       if (!axis_ref.isEmpty())
1334       {
1335         std::vector<StdString> axisList = cgPtr->getAxisList();
1336         it = std::find(axisList.begin(), axisList.end(), axis_ref.getValue());
1337         if (axisList.end() != it) domAxisScalarIds_[1] = *it;
1338       }
1339
1340       if (!scalar_ref.isEmpty())
1341       {
1342         std::vector<StdString> scalarList = cgPtr->getScalarList();
1343         it = std::find(scalarList.begin(), scalarList.end(), scalar_ref.getValue());
1344         if (scalarList.end() != it) domAxisScalarIds_[2] = *it;
1345       }
1346     }
1347     return (domAxisScalarIds_);
1348   }
1349
1350   CVariable* CField::addVariable(const string& id)
1351   {
1352     return vVariableGroup->createChild(id);
1353   }
1354
1355   CVariableGroup* CField::addVariableGroup(const string& id)
1356   {
1357     return vVariableGroup->createChildGroup(id);
1358   }
1359
1360   void CField::sendAddAllVariables()
1361   {
1362     std::vector<CVariable*> allVar = getAllVariables();
1363     std::vector<CVariable*>::const_iterator it = allVar.begin();
1364     std::vector<CVariable*>::const_iterator itE = allVar.end();
1365
1366     for (; it != itE; ++it)
1367     {
1368       this->sendAddVariable((*it)->getId());
1369       (*it)->sendAllAttributesToServer();
1370       (*it)->sendValue();
1371     }
1372   }
1373
1374   void CField::sendAddVariable(const string& id)
1375   {
1376    CContext* context = CContext::getCurrent();
1377
1378    if (!context->hasServer)
1379    {
1380       CContextClient* client = context->client;
1381
1382       CEventClient event(this->getType(),EVENT_ID_ADD_VARIABLE);
1383       if (client->isServerLeader())
1384       {
1385         CMessage msg;
1386         msg << this->getId();
1387         msg << id;
1388         const std::list<int>& ranks = client->getRanksServerLeader();
1389         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1390           event.push(*itRank,1,msg);
1391         client->sendEvent(event);
1392       }
1393       else client->sendEvent(event);
1394    }
1395   }
1396
1397   void CField::sendAddVariableGroup(const string& id)
1398   {
1399    CContext* context = CContext::getCurrent();
1400    if (!context->hasServer)
1401    {
1402       CContextClient* client = context->client;
1403
1404       CEventClient event(this->getType(),EVENT_ID_ADD_VARIABLE_GROUP);
1405       if (client->isServerLeader())
1406       {
1407         CMessage msg;
1408         msg << this->getId();
1409         msg << id;
1410         const std::list<int>& ranks = client->getRanksServerLeader();
1411         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1412           event.push(*itRank,1,msg);
1413         client->sendEvent(event);
1414       }
1415       else client->sendEvent(event);
1416    }
1417   }
1418
1419   void CField::recvAddVariable(CEventServer& event)
1420   {
1421
1422      CBufferIn* buffer = event.subEvents.begin()->buffer;
1423      string id;
1424      *buffer >> id;
1425      get(id)->recvAddVariable(*buffer);
1426   }
1427
1428   void CField::recvAddVariable(CBufferIn& buffer)
1429   {
1430      string id;
1431      buffer >> id;
1432      addVariable(id);
1433   }
1434
1435   void CField::recvAddVariableGroup(CEventServer& event)
1436   {
1437
1438      CBufferIn* buffer = event.subEvents.begin()->buffer;
1439      string id;
1440      *buffer >> id;
1441      get(id)->recvAddVariableGroup(*buffer);
1442   }
1443
1444   void CField::recvAddVariableGroup(CBufferIn& buffer)
1445   {
1446      string id;
1447      buffer >> id;
1448      addVariableGroup(id);
1449   }
1450
1451   /*!
1452    * Returns string arithmetic expression associated to the field.
1453    * \return if content is defined return content string, otherwise, if "expr" attribute is defined, return expr string.
1454    */
1455   const string& CField::getExpression(void)
1456   {
1457     if (!expr.isEmpty() && content.empty())
1458     {
1459       content = expr;
1460       expr.reset();
1461     }
1462
1463     return content;
1464   }
1465
1466   bool CField::hasExpression(void) const
1467   {
1468     return (!expr.isEmpty() || !content.empty());
1469   }
1470
1471   DEFINE_REF_FUNC(Field,field)
1472} // namespace xios
Note: See TracBrowser for help on using the repository browser.