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

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

branch re-merged with trunk @1156

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