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

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

bug fixed in mpi_comm_split. Key needs to be specifify.

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