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

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

branch merged with trunk r1130

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