source: XIOS/dev/branch_openmp/src/node/field.cpp @ 1342

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

log files

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