source: XIOS/trunk/src/node/field.cpp @ 1013

Last change on this file since 1013 was 1013, checked in by rlacroix, 7 years ago

Improve the prefetching of the fields read from a file.

Previously the data was requested at the beginning of the timestep where it would be used which was too late and could cause performance issue.

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