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

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

Inputs: Avoid sending requests for next record if we know EOF has already been reached.

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