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

Last change on this file since 1286 was 1286, checked in by oabramkina, 7 years ago

Backporting r1280 to trunk: adding a check on the field attributes freq_op and freq_offset and changing their default values.

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