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

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

The workflow is now triggered when using xios_recv_field for fields in read mode received from the servers.

Previously the workflow was triggered upon receiving the data which could cause deadlocks since there are no garanties that clients are receiving data at the same time.

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