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

Last change on this file since 1328 was 1328, checked in by yushan, 6 years ago

dev_omp

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