source: XIOS/trunk/src/node/field.cpp @ 1119

Last change on this file since 1119 was 1119, checked in by rlacroix, 7 years ago

Add the ability to check if an output field is active at the current timestep.

The "xios_field_is_active" function can now take an optional logical argument that can be used to enable this new behavior.

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