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

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

Add a new field attribut "check_if_active".

If this attribut is set to "true", xios_send_field will check if the field is active at current timestep and will be a no-op if it is not.

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