source: XIOS/dev/branch_yushan_merged/src/node/field.cpp @ 1155

Last change on this file since 1155 was 1155, checked in by yushan, 7 years ago

test_remap OK with openmp

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