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

Last change on this file since 1018 was 1018, checked in by mhnguyen, 7 years ago

Improving missing-value processing
If detect_missing_value is activated, then all missing value will be converted to
NaN (Not-a-number) in input of data flow then they will be reconverted to missing value on output

+) Update SourceFilter?, TemporalFilter? and SpatialTransformFilter? with new processing
+) Update all transformations with new processing

Test
+) On Curie
+) Work

  • 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: 46.8 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             bool hasMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
798             double defaultValue  = hasMissingValue ? default_value : (!default_value.isEmpty() ? default_value : 0.0);
799             std::pair<boost::shared_ptr<CFilter>, boost::shared_ptr<CFilter> > filters = CSpatialTransformFilter::buildFilterGraph(gc, gridRef, grid, hasMissingValue, defaultValue);
800
801             filter->connectOutput(filters.first, 0);
802             filter = filters.second;
803           }
804         }
805
806         instantDataFilter = filter;
807       }
808       // Check if we have a reference on another field
809       else if (!field_ref.isEmpty())
810         instantDataFilter = getFieldReference(gc);
811       // Check if the data is to be read from a file
812       else if (file && !file->mode.isEmpty() && file->mode == CFile::mode_attr::read)
813         instantDataFilter = serverSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid,
814                                                                                                     freq_offset.isEmpty() ? NoneDu : freq_offset,
815                                                                                                     true));
816       else // The data might be passed from the model
817       {
818          bool ignoreMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
819          double defaultValue  = ignoreMissingValue ? default_value : (!default_value.isEmpty() ? default_value : 0.0);
820          instantDataFilter = clientSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid, NoneDu, false,
821                                                                                                      ignoreMissingValue, defaultValue));       }
822     }
823
824     // If the field data is to be read by the client or/and written to a file
825     if (enableOutput && !storeFilter && !fileWriterFilter)
826     {
827       if (!read_access.isEmpty() && read_access)
828       {
829         storeFilter = boost::shared_ptr<CStoreFilter>(new CStoreFilter(gc, CContext::getCurrent(), grid));
830         instantDataFilter->connectOutput(storeFilter, 0);
831       }
832
833       if (file && (file->mode.isEmpty() || file->mode == CFile::mode_attr::write))
834       {
835         fileWriterFilter = boost::shared_ptr<CFileWriterFilter>(new CFileWriterFilter(gc, this));
836         getTemporalDataFilter(gc, file->output_freq)->connectOutput(fileWriterFilter, 0);
837       }
838     }
839   }
840
841   /*!
842    * Returns the filter needed to handle the field reference.
843    * This method should only be called when building the filter graph corresponding to the field.
844    *
845    * \param gc the garbage collector to use
846    * \return the output pin corresponding to the field reference
847    */
848   boost::shared_ptr<COutputPin> CField::getFieldReference(CGarbageCollector& gc)
849   {
850     if (instantDataFilter || field_ref.isEmpty())
851       ERROR("COutputPin* CField::getFieldReference(CGarbageCollector& gc)",
852             "Impossible to get the field reference for a field which has already been parsed or which does not have a field_ref.");
853
854     CField* fieldRef = CField::get(field_ref);
855     fieldRef->buildFilterGraph(gc, false);
856
857     std::pair<boost::shared_ptr<CFilter>, boost::shared_ptr<CFilter> > filters;
858     // Check if a spatial transformation is needed
859     if (grid && grid != fieldRef->grid && grid->hasTransform())
860     {       
861       bool hasMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
862       double defaultValue  = hasMissingValue ? default_value : (!default_value.isEmpty() ? default_value : 0.0);                               
863       filters = CSpatialTransformFilter::buildFilterGraph(gc, fieldRef->grid, grid, hasMissingValue, defaultValue);
864     }
865     else
866       filters.first = filters.second = boost::shared_ptr<CFilter>(new CPassThroughFilter(gc));
867
868     fieldRef->getInstantDataFilter()->connectOutput(filters.first, 0);
869
870     return filters.second;
871   }
872
873   /*!
874    * Returns the filter needed to handle a self reference in the field's expression.
875    * If the needed filter does not exist, it is created, otherwise it is reused.
876    * This method should only be called when building the filter graph corresponding
877    * to the field's expression.
878    *
879    * \param gc the garbage collector to use
880    * \return the output pin corresponding to a self reference
881    */
882   boost::shared_ptr<COutputPin> CField::getSelfReference(CGarbageCollector& gc)
883   {
884     if (instantDataFilter || !hasExpression())
885       ERROR("COutputPin* CField::getSelfReference(CGarbageCollector& gc)",
886             "Impossible to add a self reference to a field which has already been parsed or which does not have an expression.");
887
888     if (!selfReferenceFilter)
889     {
890       if (file && !file->mode.isEmpty() && file->mode == CFile::mode_attr::read)
891       {
892         if (!serverSourceFilter)
893           serverSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid,
894                                                                                   freq_offset.isEmpty() ? NoneDu : freq_offset,
895                                                                                   true));
896
897         selfReferenceFilter = serverSourceFilter;
898       }
899       else if (!field_ref.isEmpty())
900       {
901         CField* fieldRef = CField::get(field_ref);
902         fieldRef->buildFilterGraph(gc, false); 
903         selfReferenceFilter = fieldRef->getInstantDataFilter();
904       }
905       else
906       {
907         if (!clientSourceFilter)
908         {
909           bool ignoreMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
910           double defaultValue  = ignoreMissingValue ? default_value : (!default_value.isEmpty() ? default_value : 0.0); 
911           clientSourceFilter = boost::shared_ptr<CSourceFilter>(new CSourceFilter(gc, grid, NoneDu, false,
912                                                                                   ignoreMissingValue, defaultValue));
913         }
914
915         selfReferenceFilter = clientSourceFilter;
916       }
917     }
918
919     return selfReferenceFilter;
920   }
921
922   /*!
923    * Returns the temporal filter corresponding to the field's temporal operation
924    * for the specified operation frequency. The filter is created if it does not
925    * exist, otherwise it is reused.
926    *
927    * \param gc the garbage collector to use
928    * \param outFreq the operation frequency, i.e. the frequency at which the output data will be computed
929    * \return the output pin corresponding to the requested temporal filter
930    */
931   boost::shared_ptr<COutputPin> CField::getTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)
932   {
933     std::map<CDuration, boost::shared_ptr<COutputPin> >::iterator it = temporalDataFilters.find(outFreq);
934
935     if (it == temporalDataFilters.end())
936     {
937       if (operation.isEmpty())
938         ERROR("void CField::getTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)",
939               << "An operation must be defined for field \"" << getId() << "\".");
940
941       if (freq_op.isEmpty())
942         freq_op.setValue(TimeStep);
943       if (freq_offset.isEmpty())
944         freq_offset.setValue(NoneDu);
945
946       const bool ignoreMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
947       
948       boost::shared_ptr<CTemporalFilter> temporalFilter(new CTemporalFilter(gc, operation,
949                                                                             CContext::getCurrent()->getCalendar()->getInitDate(),
950                                                                             freq_op, freq_offset, outFreq,
951                                                                             ignoreMissingValue, ignoreMissingValue ? default_value : 0.0));
952       instantDataFilter->connectOutput(temporalFilter, 0);
953
954       it = temporalDataFilters.insert(std::make_pair(outFreq, temporalFilter)).first;
955     }
956
957     return it->second;
958   }
959
960  /*!
961    * Returns the temporal filter corresponding to the field's temporal operation
962    * for the specified operation frequency.
963    *
964    * \param gc the garbage collector to use
965    * \param outFreq the operation frequency, i.e. the frequency at which the output data will be computed
966    * \return the output pin corresponding to the requested temporal filter
967    */
968   
969   boost::shared_ptr<COutputPin> CField::getSelfTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)
970   {
971     if (instantDataFilter || !hasExpression())
972       ERROR("COutputPin* CField::getSelfTemporalDataFilter(CGarbageCollector& gc)",
973             "Impossible to add a self reference to a field which has already been parsed or which does not have an expression.");
974
975     if (!selfReferenceFilter) getSelfReference(gc) ;
976
977     if (serverSourceFilter || clientSourceFilter)
978     {
979       if (operation.isEmpty())
980         ERROR("void CField::getSelfTemporalDataFilter(CGarbageCollector& gc, CDuration outFreq)",
981               << "An operation must be defined for field \"" << getId() << "\".");
982
983       if (freq_op.isEmpty()) freq_op.setValue(TimeStep);
984       if (freq_offset.isEmpty()) freq_offset.setValue(NoneDu);
985
986       const bool ignoreMissingValue = (!detect_missing_value.isEmpty() && !default_value.isEmpty() && detect_missing_value == true);
987
988       boost::shared_ptr<CTemporalFilter> temporalFilter(new CTemporalFilter(gc, operation,
989                                                                             CContext::getCurrent()->getCalendar()->getInitDate(),
990                                                                             freq_op, freq_offset, outFreq,
991                                                                             ignoreMissingValue, ignoreMissingValue ? default_value : 0.0));
992       selfReferenceFilter->connectOutput(temporalFilter, 0);
993       return temporalFilter ;
994     }
995     else if (!field_ref.isEmpty())
996     {
997       CField* fieldRef = CField::get(field_ref);
998       fieldRef->buildFilterGraph(gc, false); 
999       return fieldRef->getTemporalDataFilter(gc, outFreq) ;
1000     }
1001  }
1002
1003   //----------------------------------------------------------------
1004/*
1005   void CField::fromBinary(StdIStream& is)
1006   {
1007      SuperClass::fromBinary(is);
1008#define CLEAR_ATT(name_)\
1009      SuperClassAttribute::operator[](#name_)->reset()
1010
1011         CLEAR_ATT(domain_ref);
1012         CLEAR_ATT(axis_ref);
1013#undef CLEAR_ATT
1014
1015   }
1016*/
1017   //----------------------------------------------------------------
1018
1019   void CField::solveGridReference(void)
1020   {
1021      if (grid_ref.isEmpty() && domain_ref.isEmpty() && axis_ref.isEmpty() && scalar_ref.isEmpty())
1022      {
1023        ERROR("CField::solveGridReference(void)",
1024              << "A grid must be defined for field '" << getFieldOutputName() << "' .");
1025      }
1026      else if (!grid_ref.isEmpty() && (!domain_ref.isEmpty() || !axis_ref.isEmpty() || !scalar_ref.isEmpty()))
1027      {
1028        ERROR("CField::solveGridReference(void)",
1029              << "Field '" << getFieldOutputName() << "' has both a grid and a domain/axis/scalar." << std::endl
1030              << "Please define either 'grid_ref' or 'domain_ref'/'axis_ref'/'scalar_ref'.");
1031      }
1032
1033      if (grid_ref.isEmpty())
1034      {
1035        std::vector<CDomain*> vecDom;
1036        std::vector<CAxis*> vecAxis;
1037        std::vector<CScalar*> vecScalar;
1038        std::vector<int> axisDomainOrderTmp;
1039       
1040        if (!domain_ref.isEmpty())
1041        {
1042          StdString tmp = domain_ref.getValue();
1043          if (CDomain::has(domain_ref))
1044          {
1045            vecDom.push_back(CDomain::get(domain_ref));
1046            axisDomainOrderTmp.push_back(2);
1047          }
1048          else
1049            ERROR("CField::solveGridReference(void)",
1050                  << "Invalid reference to domain '" << domain_ref.getValue() << "'.");
1051        }
1052
1053        if (!axis_ref.isEmpty())
1054        {
1055          if (CAxis::has(axis_ref))
1056          {
1057            vecAxis.push_back(CAxis::get(axis_ref));
1058            axisDomainOrderTmp.push_back(1);
1059          }
1060          else
1061            ERROR("CField::solveGridReference(void)",
1062                  << "Invalid reference to axis '" << axis_ref.getValue() << "'.");
1063        }
1064
1065        if (!scalar_ref.isEmpty())
1066        {
1067          if (CScalar::has(scalar_ref))
1068          {
1069            vecScalar.push_back(CScalar::get(scalar_ref));
1070            axisDomainOrderTmp.push_back(0);
1071          }
1072          else
1073            ERROR("CField::solveGridReference(void)",
1074                  << "Invalid reference to scalar '" << scalar_ref.getValue() << "'.");
1075        }
1076       
1077        CArray<int,1> axisDomainOrder(axisDomainOrderTmp.size());
1078        for (int idx = 0; idx < axisDomainOrderTmp.size(); ++idx)
1079        {
1080          axisDomainOrder(idx) = axisDomainOrderTmp[idx];
1081        }
1082
1083        // Warning: the gridId shouldn't be set as the grid_ref since it could be inherited
1084        StdString gridId = CGrid::generateId(vecDom, vecAxis, vecScalar,axisDomainOrder);
1085        if (CGrid::has(gridId))
1086          this->grid = CGrid::get(gridId);
1087        else
1088          this->grid = CGrid::createGrid(gridId, vecDom, vecAxis, vecScalar,axisDomainOrder);
1089      }
1090      else
1091      {
1092        if (CGrid::has(grid_ref))
1093          this->grid = CGrid::get(grid_ref);
1094        else
1095          ERROR("CField::solveGridReference(void)",
1096                << "Invalid reference to grid '" << grid_ref.getValue() << "'.");
1097      }
1098   }
1099
1100   void CField::solveGridDomainAxisRef(bool checkAtt)
1101   {
1102     grid->solveDomainAxisRef(checkAtt);
1103   }
1104
1105   void CField::solveCheckMaskIndex(bool doSendingIndex)
1106   {
1107     grid->checkMaskIndex(doSendingIndex);
1108   }
1109
1110   void CField::solveTransformedGrid()
1111   {
1112     if (grid && !grid->isTransformed() && hasDirectFieldReference() && grid != getDirectFieldReference()->grid)
1113     {
1114       std::vector<CGrid*> grids;
1115       // Source grid
1116       grids.push_back(getDirectFieldReference()->grid);
1117       // Intermediate grids
1118       if (!grid_path.isEmpty())
1119       {
1120         std::string gridId;
1121         size_t start = 0, end;
1122
1123         do
1124         {
1125           end = grid_path.getValue().find(',', start);
1126           if (end != std::string::npos)
1127           {
1128             gridId = grid_path.getValue().substr(start, end - start);
1129             start = end + 1;
1130           }
1131           else
1132             gridId = grid_path.getValue().substr(start);
1133
1134           if (!CGrid::has(gridId))
1135             ERROR("void CField::solveTransformedGrid()",
1136                   << "Invalid grid_path, the grid '" << gridId << "' does not exist.");
1137
1138           grids.push_back(CGrid::get(gridId));
1139         }
1140         while (end != std::string::npos);
1141       }
1142       // Destination grid
1143       grids.push_back(grid);
1144
1145       for (size_t i = 0, count = grids.size() - 1; i < count; ++i)
1146       {
1147         CGrid *gridSrc  = grids[i];
1148         CGrid *gridDest = grids[i + 1];
1149         if (!gridDest->isTransformed())
1150           gridDest->transformGrid(gridSrc);
1151       }
1152     }
1153     else if (grid && grid->hasTransform() && !grid->isTransformed())
1154     {
1155       // Temporarily deactivate the self-transformation of grid
1156       //grid->transformGrid(grid);
1157     }
1158   }
1159
1160   void CField::solveGenerateGrid()
1161   {
1162     if (grid && !grid->isTransformed() && hasDirectFieldReference() && grid != getDirectFieldReference()->grid)
1163       grid->completeGrid(getDirectFieldReference()->grid);
1164     else
1165       grid->completeGrid();
1166   }
1167
1168   void CField::solveGridDomainAxisBaseRef()
1169   {
1170     grid->solveDomainAxisRef(false);
1171     grid->solveDomainAxisBaseRef();
1172   }
1173
1174   ///-------------------------------------------------------------------
1175
1176   template <>
1177   void CGroupTemplate<CField, CFieldGroup, CFieldAttributes>::solveRefInheritance(void)
1178   {
1179      if (this->group_ref.isEmpty()) return;
1180      StdString gref = this->group_ref.getValue();
1181
1182      if (!CFieldGroup::has(gref))
1183         ERROR("CGroupTemplate<CField, CFieldGroup, CFieldAttributes>::solveRefInheritance(void)",
1184               << "[ gref = " << gref << "]"
1185               << " invalid group name !");
1186
1187      CFieldGroup* group = CFieldGroup::get(gref);
1188      CFieldGroup* owner = CFieldGroup::get(boost::polymorphic_downcast<CFieldGroup*>(this));
1189
1190      std::vector<CField*> allChildren  = group->getAllChildren();
1191      std::vector<CField*>::iterator it = allChildren.begin(), end = allChildren.end();
1192
1193      for (; it != end; it++)
1194      {
1195         CField* child = *it;
1196         if (child->hasId()) owner->createChild()->field_ref.setValue(child->getId());
1197
1198      }
1199   }
1200
1201   void CField::scaleFactorAddOffset(double scaleFactor, double addOffset)
1202   {
1203     map<int, CArray<double,1> >::iterator it;
1204     for (it = data_srv.begin(); it != data_srv.end(); it++) it->second = (it->second - addOffset) / scaleFactor;
1205   }
1206
1207   void CField::invertScaleFactorAddOffset(double scaleFactor, double addOffset)
1208   {
1209     map<int, CArray<double,1> >::iterator it;
1210     for (it = data_srv.begin(); it != data_srv.end(); it++) it->second = it->second * scaleFactor + addOffset;
1211   }
1212
1213   void CField::outputField(CArray<double,3>& fieldOut)
1214   {
1215      map<int, CArray<double,1> >::iterator it;
1216      for (it = data_srv.begin(); it != data_srv.end(); it++)
1217      {
1218        grid->outputField(it->first, it->second, fieldOut.dataFirst());
1219      }
1220   }
1221
1222   void CField::outputField(CArray<double,2>& fieldOut)
1223   {
1224      map<int, CArray<double,1> >::iterator it;
1225      for(it=data_srv.begin();it!=data_srv.end();it++)
1226      {
1227         grid->outputField(it->first, it->second, fieldOut.dataFirst());
1228      }
1229   }
1230
1231   void CField::outputField(CArray<double,1>& fieldOut)
1232   {
1233      map<int, CArray<double,1> >::iterator it;
1234
1235      for (it = data_srv.begin(); it != data_srv.end(); it++)
1236      {
1237         grid->outputField(it->first, it->second, fieldOut.dataFirst());
1238      }
1239   }
1240
1241   void CField::inputField(CArray<double,3>& fieldOut)
1242   {
1243      map<int, CArray<double,1> >::iterator it;
1244      for (it = data_srv.begin(); it != data_srv.end(); it++)
1245      {
1246        grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1247      }
1248   }
1249
1250   void CField::inputField(CArray<double,2>& fieldOut)
1251   {
1252      map<int, CArray<double,1> >::iterator it;
1253      for(it = data_srv.begin(); it != data_srv.end(); it++)
1254      {
1255         grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1256      }
1257   }
1258
1259   void CField::inputField(CArray<double,1>& fieldOut)
1260   {
1261      map<int, CArray<double,1> >::iterator it;
1262      for (it = data_srv.begin(); it != data_srv.end(); it++)
1263      {
1264         grid->inputField(it->first, fieldOut.dataFirst(), it->second);
1265      }
1266   }
1267
1268   void CField::outputCompressedField(CArray<double,1>& fieldOut)
1269   {
1270      map<int, CArray<double,1> >::iterator it;
1271
1272      for (it = data_srv.begin(); it != data_srv.end(); it++)
1273      {
1274         grid->outputCompressedField(it->first, it->second, fieldOut.dataFirst());
1275      }
1276   }
1277
1278   ///-------------------------------------------------------------------
1279
1280   void CField::parse(xml::CXMLNode& node)
1281   {
1282      SuperClass::parse(node);
1283      if (!node.getContent(this->content))
1284      {
1285        if (node.goToChildElement())
1286        {
1287          do
1288          {
1289            if (node.getElementName() == "variable" || node.getElementName() == "variable_group") this->getVirtualVariableGroup()->parseChild(node);
1290          } while (node.goToNextElement());
1291          node.goToParentElement();
1292        }
1293      }
1294    }
1295
1296   /*!
1297     This function retrieves Id of corresponding domain_ref and axis_ref (if any)
1298   of a field. In some cases, only domain exists but axis doesn't
1299   \return pair of Domain and Axis id
1300   */
1301   const std::vector<StdString>& CField::getRefDomainAxisIds()
1302   {
1303     CGrid* cgPtr = getRelGrid();
1304     if (NULL != cgPtr)
1305     {
1306       std::vector<StdString>::iterator it;
1307       if (!domain_ref.isEmpty())
1308       {
1309         std::vector<StdString> domainList = cgPtr->getDomainList();
1310         it = std::find(domainList.begin(), domainList.end(), domain_ref.getValue());
1311         if (domainList.end() != it) domAxisScalarIds_[0] = *it;
1312       }
1313
1314       if (!axis_ref.isEmpty())
1315       {
1316         std::vector<StdString> axisList = cgPtr->getAxisList();
1317         it = std::find(axisList.begin(), axisList.end(), axis_ref.getValue());
1318         if (axisList.end() != it) domAxisScalarIds_[1] = *it;
1319       }
1320
1321       if (!scalar_ref.isEmpty())
1322       {
1323         std::vector<StdString> scalarList = cgPtr->getScalarList();
1324         it = std::find(scalarList.begin(), scalarList.end(), scalar_ref.getValue());
1325         if (scalarList.end() != it) domAxisScalarIds_[2] = *it;
1326       }
1327     }
1328     return (domAxisScalarIds_);
1329   }
1330
1331   CVariable* CField::addVariable(const string& id)
1332   {
1333     return vVariableGroup->createChild(id);
1334   }
1335
1336   CVariableGroup* CField::addVariableGroup(const string& id)
1337   {
1338     return vVariableGroup->createChildGroup(id);
1339   }
1340
1341   void CField::sendAddAllVariables()
1342   {
1343     std::vector<CVariable*> allVar = getAllVariables();
1344     std::vector<CVariable*>::const_iterator it = allVar.begin();
1345     std::vector<CVariable*>::const_iterator itE = allVar.end();
1346
1347     for (; it != itE; ++it)
1348     {
1349       this->sendAddVariable((*it)->getId());
1350       (*it)->sendAllAttributesToServer();
1351       (*it)->sendValue();
1352     }
1353   }
1354
1355   void CField::sendAddVariable(const string& id)
1356   {
1357    CContext* context = CContext::getCurrent();
1358
1359    if (!context->hasServer)
1360    {
1361       CContextClient* client = context->client;
1362
1363       CEventClient event(this->getType(),EVENT_ID_ADD_VARIABLE);
1364       if (client->isServerLeader())
1365       {
1366         CMessage msg;
1367         msg << this->getId();
1368         msg << id;
1369         const std::list<int>& ranks = client->getRanksServerLeader();
1370         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1371           event.push(*itRank,1,msg);
1372         client->sendEvent(event);
1373       }
1374       else client->sendEvent(event);
1375    }
1376   }
1377
1378   void CField::sendAddVariableGroup(const string& id)
1379   {
1380    CContext* context = CContext::getCurrent();
1381    if (!context->hasServer)
1382    {
1383       CContextClient* client = context->client;
1384
1385       CEventClient event(this->getType(),EVENT_ID_ADD_VARIABLE_GROUP);
1386       if (client->isServerLeader())
1387       {
1388         CMessage msg;
1389         msg << this->getId();
1390         msg << id;
1391         const std::list<int>& ranks = client->getRanksServerLeader();
1392         for (std::list<int>::const_iterator itRank = ranks.begin(), itRankEnd = ranks.end(); itRank != itRankEnd; ++itRank)
1393           event.push(*itRank,1,msg);
1394         client->sendEvent(event);
1395       }
1396       else client->sendEvent(event);
1397    }
1398   }
1399
1400   void CField::recvAddVariable(CEventServer& event)
1401   {
1402
1403      CBufferIn* buffer = event.subEvents.begin()->buffer;
1404      string id;
1405      *buffer >> id;
1406      get(id)->recvAddVariable(*buffer);
1407   }
1408
1409   void CField::recvAddVariable(CBufferIn& buffer)
1410   {
1411      string id;
1412      buffer >> id;
1413      addVariable(id);
1414   }
1415
1416   void CField::recvAddVariableGroup(CEventServer& event)
1417   {
1418
1419      CBufferIn* buffer = event.subEvents.begin()->buffer;
1420      string id;
1421      *buffer >> id;
1422      get(id)->recvAddVariableGroup(*buffer);
1423   }
1424
1425   void CField::recvAddVariableGroup(CBufferIn& buffer)
1426   {
1427      string id;
1428      buffer >> id;
1429      addVariableGroup(id);
1430   }
1431
1432   /*!
1433    * Returns string arithmetic expression associated to the field.
1434    * \return if content is defined return content string, otherwise, if "expr" attribute is defined, return expr string.
1435    */
1436   const string& CField::getExpression(void)
1437   {
1438     if (!expr.isEmpty() && content.empty())
1439     {
1440       content = expr;
1441       expr.reset();
1442     }
1443
1444     return content;
1445   }
1446
1447   bool CField::hasExpression(void) const
1448   {
1449     return (!expr.isEmpty() || !content.empty());
1450   }
1451
1452   DEFINE_REF_FUNC(Field,field)
1453} // namespace xios
Note: See TracBrowser for help on using the repository browser.