Salome HOME
Improve the clearness of the message in case of the ELNO surface filter is used inste...
[modules/paravis.git] / src / Plugins / MEDReader / IO / MEDFileFieldRepresentationTree.cxx
1 // Copyright (C) 2010-2014  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 // Author : Anthony Geay
20
21 #include "MEDTimeReq.hxx"
22 #include "MEDUtilities.hxx"
23
24 #include "MEDFileFieldRepresentationTree.hxx"
25 #include "MEDCouplingFieldDiscretization.hxx"
26 #include "MEDCouplingFieldDouble.hxx"
27 #include "InterpKernelGaussCoords.hxx"
28 #include "MEDFileData.hxx"
29 #include "SauvReader.hxx"
30
31 #include "vtkXMLUnstructuredGridWriter.h"//
32
33 #include "vtkUnstructuredGrid.h"
34 #include "vtkRectilinearGrid.h"
35 #include "vtkStructuredGrid.h"
36 #include "vtkUnsignedCharArray.h"
37 #include "vtkQuadratureSchemeDefinition.h"
38 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
39 #include "vtkInformationIntegerKey.h"
40 #include "vtkInformation.h"
41 #include "vtkIdTypeArray.h"
42 #include "vtkDoubleArray.h"
43 #include "vtkIntArray.h"
44 #include "vtkCellArray.h"
45 #include "vtkPointData.h"
46 #include "vtkFieldData.h"
47 #include "vtkCellData.h"
48
49 #include "vtksys/stl/string"
50 #include "vtksys/ios/fstream"
51 #include "vtksys/stl/algorithm"
52 #include "vtkMutableDirectedGraph.h"
53
54 using namespace ParaMEDMEM;
55
56 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
57
58 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
59
60 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
61
62 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
63
64 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
65
66 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
67
68 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
69
70 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
71
72 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
73
74 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
75
76 vtkIdTypeArray *ELGACmp::findOrCreate(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, bool& isNew) const
77 {
78   vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
79   if(try0)
80     {
81       isNew=false;
82       return try0;
83     }
84   else
85     {
86       isNew=true;
87       return createNew(globs,locsReallyUsed,vtkd,ds);
88     }
89 }
90
91 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd) const
92 {
93   std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
94   if(it==_loc_names.end())
95     return 0;
96   std::size_t pos(std::distance(_loc_names.begin(),it));
97   vtkIdTypeArray *ret(_elgas[pos]);
98   vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
99   for(std::vector<std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it=_defs[pos].begin();it!=_defs[pos].end();it++)
100     {
101       key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
102     }
103   vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
104   return ret;
105 }
106
107 vtkIdTypeArray *ELGACmp::createNew(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds) const
108 {
109   static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
110   std::vector< std::vector<std::string> > locNames(_loc_names);
111   std::vector<vtkIdTypeArray *> elgas(_elgas);
112   std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > defs;
113   //
114   std::vector< std::vector<std::string> >::const_iterator it(std::find(locNames.begin(),locNames.end(),locsReallyUsed));
115   if(it!=locNames.end())
116     throw INTERP_KERNEL::Exception("ELGACmp::createNew : Method is expected to be called after isExisting call ! Entry already exists !");
117   locNames.push_back(locsReallyUsed);
118   vtkIdTypeArray *elga(vtkIdTypeArray::New());
119   elga->SetNumberOfComponents(1);
120   vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
121   std::map<unsigned char,int> m;
122   for(std::vector<std::string>::const_iterator it=locsReallyUsed.begin();it!=locsReallyUsed.end();it++)
123     {
124       vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
125       const MEDFileFieldLoc& loc(globs->getLocalization((*it).c_str()));
126       INTERP_KERNEL::NormalizedCellType ct(loc.getGeoType());
127       const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
128       int nbGaussPt(loc.getNbOfGaussPtPerCell()),nbPtsPerCell((int)cm.getNumberOfNodes()),dimLoc(loc.getDimension());
129       // WARNING : these 2 lines are a workaround, due to users that write a ref element with dimension not equal to dimension of the geometric element.
130       std::vector<double> gsCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getGaussCoords()));
131       std::vector<double> refCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getRefCoords()));
132       double *shape(new double[nbPtsPerCell*nbGaussPt]);
133       INTERP_KERNEL::GaussInfo calculator(ct,gsCoods2,nbGaussPt,refCoods2,nbPtsPerCell);
134       calculator.initLocalInfo();
135       const std::vector<double>& wgths(loc.getGaussWeights());
136       for(int i=0;i<nbGaussPt;i++)
137         {
138           const double *pt0(calculator.getFunctionValues(i));
139           std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
140         }
141       unsigned char vtkType(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE[ct]);
142       m[vtkType]=nbGaussPt;
143       def->Initialize(vtkType,nbPtsPerCell,nbGaussPt,shape,const_cast<double *>(&wgths[0]));
144       delete [] shape;
145       key->Set(elga->GetInformation(),def,vtkType);
146       key->Set(vtkd->GetInformation(),def,vtkType);
147       defs.push_back(std::pair< vtkQuadratureSchemeDefinition *, unsigned char >(def,vtkType));
148     }
149   //
150   vtkIdType ncell(ds->GetNumberOfCells());
151   int *pt(new int[ncell]),offset(0);
152   for(vtkIdType cellId=0;cellId<ncell;cellId++)
153     {
154       vtkCell *cell(ds->GetCell(cellId));
155       int delta(m[cell->GetCellType()]);
156       pt[cellId]=offset;
157       offset+=delta;
158     }
159   elga->GetInformation()->Set(MEDUtilities::ELGA(),1);
160   elga->SetArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
161   std::ostringstream oss; oss << "ELGA" << "@" << _loc_names.size();
162   std::string ossStr(oss.str());
163   elga->SetName(ossStr.c_str());
164   elga->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
165   vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elga->GetName());
166   elgas.push_back(elga);
167   //
168   _loc_names=locNames;
169   _elgas=elgas;
170   _defs.push_back(defs);
171 }
172
173 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
174 {
175   for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
176     ds->GetCellData()->AddArray(*it);
177 }
178
179 ELGACmp::~ELGACmp()
180 {
181   for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
182     (*it)->Delete();
183   for(std::vector< std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > >::const_iterator it0=_defs.begin();it0!=_defs.end();it0++)
184     for(std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
185       (*it1).first->Delete();
186 }
187
188 //=
189
190 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
191 {
192 }
193
194 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>& arr):ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
195 {
196   std::vector< std::vector<ParaMEDMEM::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
197   if(typs.size()<1)
198     throw INTERP_KERNEL::Exception("There is a big internal problem in MEDLoader ! The field time spitting has failed ! A CRASH will occur soon !");
199   if(typs[0].size()!=1)
200     throw INTERP_KERNEL::Exception("There is a big internal problem in MEDLoader ! The field spitting by spatial discretization has failed ! A CRASH will occur soon !");
201   ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDiscretization> fd(MEDCouplingFieldDiscretization::New(typs[0][0]));
202   std::ostringstream oss2; oss2 << (operator->())->getName() << ZE_SEP << fd->getRepr();
203   _ze_name=oss2.str();
204 }
205
206 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
207 {
208   ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>::operator=(other);
209   _id=-1;
210   _activated=false;
211   _ze_name=other._ze_name;
212   _ze_full_name.clear();
213   return *this;
214 }
215
216 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
217 {
218   _id=id++;
219 }
220
221 int MEDFileFieldRepresentationLeavesArrays::getId() const
222 {
223   return _id;
224 }
225
226 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
227 {
228   return _ze_full_name;
229 }
230
231 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, const std::string& tsName, const std::string& meshName, const std::string& comSupStr, std::vector<std::string>& names) const
232 {
233   vtkIdType refId(sil->AddChild(root,edge));
234   names.push_back(_ze_name);
235   std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
236   _ze_full_name=oss3.str();
237   //
238   if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
239     {
240       sil->AddChild(refId,edge);
241       names.push_back(std::string());
242     }
243 }
244
245 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
246 {
247   return _activated;
248 }
249
250 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
251 {
252   bool ret(_activated!=status);
253   _activated=status;
254   return ret;
255 }
256
257 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshStruct *mst, vtkDataSet *ds) const
258 {
259   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
260   static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
261   tr->setNumberOfTS((operator->())->getNumberOfTS());
262   tr->initIterator();
263   for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
264     {
265       MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeField1TS> f1ts((operator->())->getTimeStepAtPos(tr->getCurrent()));
266       MEDFileAnyTypeField1TS *f1tsPtr(f1ts);
267       MEDFileField1TS *f1tsPtrDbl(dynamic_cast<MEDFileField1TS *>(f1tsPtr));
268       MEDFileIntField1TS *f1tsPtrInt(dynamic_cast<MEDFileIntField1TS *>(f1tsPtr));
269       DataArray *crudeArr(0),*postProcessedArr(0);
270       if(f1tsPtrDbl)
271         crudeArr=f1tsPtrDbl->getUndergroundDataArray();
272       else if(f1tsPtrInt)
273         crudeArr=f1tsPtrInt->getUndergroundDataArray();
274       else
275         throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment !");
276       MEDFileField1TSStructItem fsst(MEDFileField1TSStructItem::BuildItemFrom(f1ts,mst));
277       f1ts->loadArraysIfNecessary();
278       MEDCouplingAutoRefCountObjectPtr<DataArray> v(mml->buildDataArray(fsst,globs,crudeArr));
279       postProcessedArr=v;
280       //
281       std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
282       if(discs.size()!=1)
283         throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
284       vtkFieldData *att(0);
285       switch(discs[0])
286         {
287         case ON_CELLS:
288           {
289             att=ds->GetCellData();
290             break;
291           }
292         case ON_NODES:
293           {
294             att=ds->GetPointData();
295             break;
296           }
297         case ON_GAUSS_NE:
298           {
299             att=ds->GetFieldData();
300             break;
301           }
302         case ON_GAUSS_PT:
303           {
304             att=ds->GetFieldData();
305             break;
306           }
307         default:
308           throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
309         }
310       if(f1tsPtrDbl)
311         {
312           DataArray *vPtr(v); DataArrayDouble *vd(static_cast<DataArrayDouble *>(vPtr));
313           vtkDoubleArray *vtkd(vtkDoubleArray::New());
314           vtkd->SetNumberOfComponents(vd->getNumberOfComponents());
315           for(int i=0;i<vd->getNumberOfComponents();i++)
316             vtkd->SetComponentName(i,vd->getInfoOnComponent(i).c_str());
317           if(postProcessedArr!=crudeArr)
318             {
319               vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vd->accessToMemArray().setSpecificDeallocator(0);
320             }
321           else
322             {
323               vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
324             }
325           std::string name(tr->buildName(f1ts->getName()));
326           vtkd->SetName(name.c_str());
327           att->AddArray(vtkd);
328           vtkd->Delete();
329           if(discs[0]==ON_GAUSS_PT)
330             {
331               bool tmp;
332               _elga_cmp.findOrCreate(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp);
333             }
334           if(discs[0]==ON_GAUSS_NE)
335             {
336               vtkIdTypeArray *elno(vtkIdTypeArray::New());
337               elno->SetNumberOfComponents(1);
338               vtkIdType ncell(ds->GetNumberOfCells());
339               int *pt(new int[ncell]),offset(0);
340               std::set<int> cellTypes;
341               for(vtkIdType cellId=0;cellId<ncell;cellId++)
342                 {
343                   vtkCell *cell(ds->GetCell(cellId));
344                   int delta(cell->GetNumberOfPoints());
345                   cellTypes.insert(cell->GetCellType());
346                   pt[cellId]=offset;
347                   offset+=delta;
348                 }
349               elno->GetInformation()->Set(MEDUtilities::ELNO(),1);
350               elno->SetArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
351               std::string nameElno("ELNO"); nameElno+="@"; nameElno+=name;
352               elno->SetName(nameElno.c_str());
353               ds->GetCellData()->AddArray(elno);
354               vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elno->GetName());
355               elno->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
356               //
357               vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
358               for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
359                 {
360                   const unsigned char *pos(std::find(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,*it));
361                   if(pos==MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
362                     continue;
363                   INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)std::distance(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos));
364                   const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
365                   int nbGaussPt(cm.getNumberOfNodes()),dim(cm.getDimension());
366                   vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
367                   double *shape(new double[nbGaussPt*nbGaussPt]);
368                   std::size_t dummy;
369                   const double *gsCoords(MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType(ct,dummy));
370                   const double *refCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));
371                   const double *weights(MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(ct,dummy));
372                   std::vector<double> gsCoords2(gsCoords,gsCoords+nbGaussPt*dim),refCoords2(refCoords,refCoords+nbGaussPt*dim);
373                   INTERP_KERNEL::GaussInfo calculator(ct,gsCoords2,nbGaussPt,refCoords2,nbGaussPt);
374                   calculator.initLocalInfo();
375                   for(int i=0;i<nbGaussPt;i++)
376                     {
377                       const double *pt0(calculator.getFunctionValues(i));
378                       std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
379                     }
380                   def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
381                   delete [] shape;
382                   key->Set(elno->GetInformation(),def,*it);
383                   key->Set(vtkd->GetInformation(),def,*it);
384                   def->Delete();
385                 }
386               //
387               elno->Delete();
388             }
389         }
390       else if(f1tsPtrInt)
391         {
392           DataArray *vPtr(v); DataArrayInt *vi(static_cast<DataArrayInt *>(vPtr));
393           vtkIntArray *vtkd(vtkIntArray::New());
394           vtkd->SetNumberOfComponents(vi->getNumberOfComponents());
395           for(int i=0;i<vi->getNumberOfComponents();i++)
396             vtkd->SetComponentName(i,vi->getVarOnComponent(i).c_str());
397           if(postProcessedArr!=crudeArr)
398             {
399               vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vi->accessToMemArray().setSpecificDeallocator(0);
400             }
401           else
402             {
403               vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
404             }
405           std::string name(tr->buildName(f1ts->getName()));
406           vtkd->SetName(name.c_str());
407           att->AddArray(vtkd);
408           vtkd->Delete();
409         }
410       else
411         throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
412     }
413 }
414
415 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
416 {
417   _elga_cmp.appendELGAIfAny(ds);
418 }
419
420 ////////////////////
421
422 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
423 {
424 }
425
426 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves(const std::vector< ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS> >& arr,
427                                                                    const ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFastCellSupportComparator>& fsp):_arrays(arr.size()),_fsp(fsp),_cached_ds(0)
428 {
429   for(std::size_t i=0;i<arr.size();i++)
430     _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
431 }
432
433 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
434 {
435   if(_cached_ds)
436     _cached_ds->Delete();
437 }
438
439 bool MEDFileFieldRepresentationLeaves::empty() const
440 {
441   const MEDFileFastCellSupportComparator *fcscp(_fsp);
442   return fcscp==0 || _arrays.empty();
443 }
444
445 void MEDFileFieldRepresentationLeaves::setId(int& id) const
446 {
447   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
448     (*it).setId(id);
449 }
450
451 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
452 {
453   return _arrays[0]->getMeshName();
454 }
455
456 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
457 {
458   return (int)_arrays.size();
459 }
460
461 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
462 {
463   return _arrays[0]->getNumberOfTS();
464 }
465
466 /*!
467  * \param [in] ms is the meshes pointer. It can be used only for information of geometric types. No special processing will be requested on ms.
468  */
469 void MEDFileFieldRepresentationLeaves::feedSIL(const ParaMEDMEM::MEDFileMeshes *ms, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, const std::string& tsName, const std::string& meshName, const std::string& comSupStr, std::vector<std::string>& names) const
470 {
471   vtkIdType root2(sil->AddChild(root,edge));
472   names.push_back(std::string("Arrs"));
473   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
474     (*it).feedSIL(sil,root2,edge,tsName,meshName,comSupStr,names);
475   //
476   vtkIdType root3(sil->AddChild(root,edge));
477   names.push_back(std::string("InfoOnGeoType"));
478   const ParaMEDMEM::MEDFileMesh *m(0);
479   if(ms)
480     m=ms->getMeshWithName(meshName);
481   const ParaMEDMEM::MEDFileFastCellSupportComparator *fsp(_fsp);
482   if(!fsp || fsp->getNumberOfTS()==0)
483     return ;
484   std::vector< INTERP_KERNEL::NormalizedCellType > gts(fsp->getGeoTypesAt(0,m));
485   for(std::vector< INTERP_KERNEL::NormalizedCellType >::const_iterator it2=gts.begin();it2!=gts.end();it2++)
486     {
487       const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(*it2));
488       std::string cmStr(cm.getRepr()); cmStr=cmStr.substr(5);//skip "NORM_"
489       sil->AddChild(root3,edge);
490       names.push_back(cmStr);
491     }
492 }
493
494 bool MEDFileFieldRepresentationLeaves::containId(int id) const
495 {
496   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
497     if((*it).getId()==id)
498       return true;
499   return false;
500 }
501
502 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
503 {
504   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
505     if((*it).getZeName()==name)
506       {
507         id=(*it).getId();
508         return true;
509       }
510   return false;
511 }
512
513 bool MEDFileFieldRepresentationLeaves::isActivated() const
514 {
515   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
516     if((*it).getStatus())
517       return true;
518   return false;
519 }
520
521 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
522 {
523   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
524     {
525       os << "         - " << (*it0).getZeName() << " (";
526       if((*it0).getStatus())
527         os << "X";
528       else
529         os << " ";
530       os << ")" << std::endl;
531     }
532 }
533
534 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
535 {
536   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
537     (*it).setStatus(true);
538 }
539
540 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
541 {
542   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
543     if((*it).getId()==id)
544       return *it;
545   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
546 }
547
548 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
549 {
550   if(_arrays.size()<1)
551     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
552   std::vector<double> ret;
553   std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
554   return tk.getTimeStepsRegardingPolicy(dtits,ret);
555 }
556
557 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
558 {
559   if(!_arrays.empty())
560     return _arrays[0]->getTimeSteps(ts);
561   else
562     {
563       ts.clear();
564       return std::vector< std::pair<int,int> >();
565     }
566 }
567
568 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
569 {
570   std::ostringstream oss;
571   oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
572   std::vector<double> ret1;
573   std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
574   std::size_t sz(ret1.size());
575   for(std::size_t i=0;i<sz;i++)
576     {
577       oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
578       if(i!=sz-1)
579         oss << ", ";
580       std::string tmp(oss.str());
581       if(tmp.size()>200 && i!=sz-1)
582         {
583           oss << "...";
584           break;
585         }
586     }
587   oss << ")";
588   return oss.str();
589 }
590
591 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshes *meshes, vtkDataSet *ds) const
592 {
593   if(_arrays.size()<1)
594     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
595   MEDCouplingAutoRefCountObjectPtr<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
596   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
597     if((*it).getStatus())
598       {
599         (*it).appendFields(tr,globs,mml,mst,ds);
600         (*it).appendELGAIfAny(ds);
601       }
602 }
603
604 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
605 {
606   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
607   DataArrayDouble *coordsMC(0);
608   DataArrayByte *typesMC(0);
609   DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
610   bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
611   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coordsSafe(coordsMC);
612   MEDCouplingAutoRefCountObjectPtr<DataArrayByte> typesSafe(typesMC);
613   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
614   //
615   int nbOfCells(typesSafe->getNbOfElems());
616   vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
617   vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
618   cellTypes->SetArray(reinterpret_cast<unsigned char *>(typesSafe->getPointer()),nbOfCells,0,VTK_DATA_ARRAY_FREE); typesSafe->accessToMemArray().setSpecificDeallocator(0);
619   vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
620   cellLocations->SetArray(cellLocationsSafe->getPointer(),nbOfCells,0,VTK_DATA_ARRAY_FREE); cellLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
621   vtkCellArray *cells(vtkCellArray::New());
622   vtkIdTypeArray *cells2(vtkIdTypeArray::New());
623   cells2->SetArray(cellsSafe->getPointer(),cellsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); cellsSafe->accessToMemArray().setSpecificDeallocator(0);
624   cells->SetCells(nbOfCells,cells2);
625   cells2->Delete();
626   if(faceLocationsMC!=0 && facesMC!=0)
627     {
628       vtkIdTypeArray *faces(vtkIdTypeArray::New());
629       faces->SetArray(facesSafe->getPointer(),facesSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); facesSafe->accessToMemArray().setSpecificDeallocator(0);
630       vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
631       faceLocations->SetArray(faceLocationsSafe->getPointer(),faceLocationsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); faceLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
632       ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
633       faceLocations->Delete();
634       faces->Delete();
635     }
636   else
637     ret->SetCells(cellTypes,cellLocations,cells);
638   cellTypes->Delete();
639   cellLocations->Delete();
640   cells->Delete();
641   vtkPoints *pts(vtkPoints::New());
642   vtkDoubleArray *pts2(vtkDoubleArray::New());
643   pts2->SetNumberOfComponents(3);
644   if(!statusOfCoords)
645     {
646       pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);
647       coordsSafe->accessToMemArray().setSpecificDeallocator(0);
648     }
649   else
650     pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
651   pts->SetData(pts2);
652   pts2->Delete();
653   ret->SetPoints(pts);
654   pts->Delete();
655   //
656   return ret;
657 }
658
659 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(ParaMEDMEM::MEDCMeshMultiLev *mm) const
660 {
661   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
662   bool isInternal;
663   std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
664   vtkDoubleArray *vtkTmp(0);
665   vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
666   std::size_t dim(arrs.size());
667   if(dim<1 || dim>3)
668     throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
669   int sizePerAxe[3]={1,1,1};
670   sizePerAxe[0]=arrs[0]->getNbOfElems();
671   if(dim>=2)
672     sizePerAxe[1]=arrs[1]->getNbOfElems();
673   if(dim==3)
674     sizePerAxe[2]=arrs[2]->getNbOfElems();
675   ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
676   vtkTmp=vtkDoubleArray::New();
677   vtkTmp->SetNumberOfComponents(1);
678   if(isInternal)
679     vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
680   else
681     { vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[0]->accessToMemArray().setSpecificDeallocator(0); }
682   ret->SetXCoordinates(vtkTmp);
683   vtkTmp->Delete();
684   arrs[0]->decrRef();
685   if(dim>=2)
686     {
687       vtkTmp=vtkDoubleArray::New();
688       vtkTmp->SetNumberOfComponents(1);
689       if(isInternal)
690         vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
691       else
692         { vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[1]->accessToMemArray().setSpecificDeallocator(0); }
693       ret->SetYCoordinates(vtkTmp);
694       vtkTmp->Delete();
695       arrs[1]->decrRef();
696     }
697   if(dim==3)
698     {
699       vtkTmp=vtkDoubleArray::New();
700       vtkTmp->SetNumberOfComponents(1);
701       if(isInternal)
702         vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
703       else
704         { vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[2]->accessToMemArray().setSpecificDeallocator(0); }
705       ret->SetZCoordinates(vtkTmp);
706       vtkTmp->Delete();
707       arrs[2]->decrRef();
708     }
709   return ret;
710 }
711
712 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(ParaMEDMEM::MEDCurveLinearMeshMultiLev *mm) const
713 {
714   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
715   int meshStr[3]={1,1,1};
716   DataArrayDouble *coords(0);
717   std::vector<int> nodeStrct;
718   bool isInternal;
719   mm->buildVTUArrays(coords,nodeStrct,isInternal);
720   std::size_t dim(nodeStrct.size());
721   if(dim<1 || dim>3)
722     throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
723   meshStr[0]=nodeStrct[0];
724   if(dim>=2)
725     meshStr[1]=nodeStrct[1];
726   if(dim==3)
727     meshStr[2]=nodeStrct[2];
728   vtkStructuredGrid *ret(vtkStructuredGrid::New());
729   ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
730   vtkDoubleArray *da(vtkDoubleArray::New());
731   da->SetNumberOfComponents(3);
732   if(coords->getNumberOfComponents()==3)
733     {
734       if(isInternal)
735         da->SetArray(coords->getPointer(),coords->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);//VTK has not the ownership of double * because MEDLoader main struct has it !
736       else
737         { da->SetArray(coords->getPointer(),coords->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); coords->accessToMemArray().setSpecificDeallocator(0); }
738     }
739   else
740     {
741       MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
742       da->SetArray(coords2->getPointer(),coords2->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);//let VTK deal with double *
743       coords2->accessToMemArray().setSpecificDeallocator(0);
744     }
745   coords->decrRef();
746   vtkPoints *points=vtkPoints::New();
747   ret->SetPoints(points);
748   points->SetData(da);
749   points->Delete();
750   da->Delete();
751   return ret;
752 }
753  
754 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDFileMeshes *meshes) const
755 {
756   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
757   vtkDataSet *ret(0);
758   //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
759   MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change 
760   MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml2(mml->prepare());
761   MEDMeshMultiLev *ptMML2(mml2);
762   if(!_cached_ds)
763     {
764       MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
765       MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
766       MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
767       
768       if(ptUMML2)
769         {
770           ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
771         }
772       else if(ptCMML2)
773         {
774           ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
775         }
776       else if(ptCLMML2)
777         {
778           ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
779         }
780       else
781         throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
782       _cached_ds=ret->NewInstance();
783       _cached_ds->ShallowCopy(ret);
784     }
785   else
786     {
787       ret=_cached_ds->NewInstance();
788       ret->ShallowCopy(_cached_ds);
789     }
790   //
791   appendFields(tr,globs,mml,meshes,ret);
792   // The arrays links to mesh
793   DataArrayInt *famCells(0),*numCells(0);
794   bool noCpyFamCells(false),noCpyNumCells(false);
795   ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
796   if(famCells)
797     {
798       vtkIntArray *vtkTab(vtkIntArray::New());
799       vtkTab->SetNumberOfComponents(1);
800       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
801       if(noCpyFamCells)
802         vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
803       else
804         { vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famCells->accessToMemArray().setSpecificDeallocator(0); }
805       ret->GetCellData()->AddArray(vtkTab);
806       vtkTab->Delete();
807       famCells->decrRef();
808     }
809   ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
810   if(numCells)
811     {
812       vtkIntArray *vtkTab(vtkIntArray::New());
813       vtkTab->SetNumberOfComponents(1);
814       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
815       if(noCpyNumCells)
816         vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
817       else
818         { vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numCells->accessToMemArray().setSpecificDeallocator(0); }
819       ret->GetCellData()->AddArray(vtkTab);
820       vtkTab->Delete();
821       numCells->decrRef();
822     }
823   // The arrays links to mesh
824   DataArrayInt *famNodes(0),*numNodes(0);
825   bool noCpyFamNodes(false),noCpyNumNodes(false);
826   ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
827   if(famNodes)
828     {
829       vtkIntArray *vtkTab(vtkIntArray::New());
830       vtkTab->SetNumberOfComponents(1);
831       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
832       if(noCpyFamNodes)
833         vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
834       else
835         { vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famNodes->accessToMemArray().setSpecificDeallocator(0); }
836       ret->GetPointData()->AddArray(vtkTab);
837       vtkTab->Delete();
838       famNodes->decrRef();
839     }
840   ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
841   if(numNodes)
842     {
843       vtkIntArray *vtkTab(vtkIntArray::New());
844       vtkTab->SetNumberOfComponents(1);
845       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
846       if(noCpyNumNodes)
847         vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
848       else
849         { vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numNodes->accessToMemArray().setSpecificDeallocator(0); }
850       ret->GetPointData()->AddArray(vtkTab);
851       vtkTab->Delete();
852       numNodes->decrRef();
853     }
854   return ret;
855 }
856
857 //////////////////////
858
859 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
860 {
861 }
862
863 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
864 {
865   int ret(0);
866   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
867     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
868       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
869         ret+=(*it2).getNumberOfArrays();
870   return ret;
871 }
872
873 void MEDFileFieldRepresentationTree::assignIds() const
874 {
875   int zeId(0);
876   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
877     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
878       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
879         (*it2).setId(zeId);
880 }
881 void MEDFileFieldRepresentationTree::activateTheFirst() const
882 {
883   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
884     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
885       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
886         {
887           (*it2).activateAllArrays();
888           return ;
889         }
890 }
891
892 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
893 {
894   std::size_t it0Cnt(0);
895   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
896     {
897       vtkIdType InfoOnTSId(sil->AddChild(root,edge));
898       names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
899       //
900       vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
901       std::vector<double> ts;
902       std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
903       std::size_t nbOfTS(dtits.size());
904       std::ostringstream oss3; oss3 << nbOfTS;
905       names.push_back(oss3.str());
906       for(std::size_t i=0;i<nbOfTS;i++)
907         {
908           std::ostringstream oss4; oss4 << dtits[i].first;
909           vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
910           names.push_back(oss4.str());
911           std::ostringstream oss5; oss5 << dtits[i].second;
912           vtkIdType ItId(sil->AddChild(DtId,edge));
913           names.push_back(oss5.str());
914           std::ostringstream oss6; oss6 << ts[i];
915           sil->AddChild(ItId,edge);
916           names.push_back(oss6.str());
917         }
918       //
919       std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
920       std::string tsName(oss.str());
921       vtkIdType typeId0(sil->AddChild(root,edge));
922       names.push_back(tsName);
923       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
924         {
925           std::string meshName((*it1)[0].getMeshName());
926           vtkIdType typeId1(sil->AddChild(typeId0,edge));
927           names.push_back(meshName);
928           std::size_t it2Cnt(0);
929           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
930             {
931               std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
932               std::string comSupStr(oss2.str());
933               vtkIdType typeId2(sil->AddChild(typeId1,edge));
934               names.push_back(comSupStr);
935               (*it2).feedSIL(_ms,sil,typeId2,edge,tsName,meshName,comSupStr,names);
936             } 
937         }
938     }
939 }
940
941 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
942 {
943   int dummy0(0),dummy1(0),dummy2(0);
944   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
945   std::string ret(leaf.getMeshName());
946   int i(0);
947   MEDFileMesh *m(0);
948   for(;i<_ms->getNumberOfMeshes();i++)
949     {
950       m=_ms->getMeshAtPos(i);
951       if(m->getName()==ret)
952         break;
953     }
954   if(i==_ms->getNumberOfMeshes())
955     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
956   vtkIdType typeId0(sil->AddChild(root,edge));
957   names.push_back(m->getName());
958   //
959   vtkIdType typeId1(sil->AddChild(typeId0,edge));
960   names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
961   std::vector<std::string> grps(m->getGroupsNames());
962   for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
963     {
964       vtkIdType typeId2(sil->AddChild(typeId1,edge));
965       names.push_back(*it0);
966       std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
967       for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
968         {
969           sil->AddChild(typeId2,edge);
970           names.push_back((*it1).c_str());
971         }
972     }
973   //
974   vtkIdType typeId11(sil->AddChild(typeId0,edge));
975   names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
976   std::vector<std::string> fams(m->getFamiliesNames());
977   for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
978     {
979       sil->AddChild(typeId11,edge);
980       int famId(m->getFamilyId((*it00).c_str()));
981       std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
982       names.push_back(oss.str());
983     }
984   return ret;
985 }
986
987 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
988 {
989   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
990     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
991       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
992         if((*it2).containId(id))
993           return (*it2).getLeafArr(id);
994   throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
995 }
996
997 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
998 {
999   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1000   return elt.getZeName();
1001 }
1002
1003 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1004 {
1005   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1006   return elt.getStatus();
1007 }
1008
1009 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1010 {
1011   int ret(-1);
1012   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1013     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1014       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1015         if((*it2).containZeName(name,ret))
1016           return ret;
1017   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name !");
1018 }
1019
1020 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1021 {
1022   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1023   bool ret(elt.setStatus(status));//to be implemented
1024   return ret;
1025 }
1026
1027 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1028 {
1029   int ret(0);
1030   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1031     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1032       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1033         ret=std::max(ret,(*it2).getNumberOfTS());
1034   return ret;
1035 }
1036
1037 /*!
1038  * 
1039  */
1040 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv)
1041 {
1042   if(isMEDOrSauv)
1043     {
1044       _ms=MEDFileMeshes::New(fileName);
1045       _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
1046     }
1047   else
1048     {
1049       MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::SauvReader> sr(ParaMEDMEM::SauvReader::New(fileName));
1050       MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileData> mfd(sr->loadInMEDFileDS());
1051       _ms=mfd->getMeshes(); _ms->incrRef();
1052       int nbMeshes(_ms->getNumberOfMeshes());
1053       for(int i=0;i<nbMeshes;i++)
1054         {
1055           ParaMEDMEM::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1056           ParaMEDMEM::MEDFileUMesh *tmp2(dynamic_cast<ParaMEDMEM::MEDFileUMesh *>(tmp));
1057           if(tmp2)
1058             tmp2->forceComputationOfParts();
1059         }
1060       _fields=mfd->getFields();
1061       if((ParaMEDMEM::MEDFileFields *)_fields)
1062         _fields->incrRef();
1063     }
1064   if(!((ParaMEDMEM::MEDFileFields *)_fields))
1065     {
1066       _fields=BuildFieldFromMeshes(_ms);
1067     }
1068   else
1069     {
1070       AppendFieldFromMeshes(_ms,_fields);
1071     }
1072   std::vector<std::string> meshNames(_ms->getMeshesNames());
1073   std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> > fields_per_mesh(meshNames.size());
1074   for(std::size_t i=0;i<meshNames.size();i++)
1075     {
1076       fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1077     }
1078   std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1079   std::size_t k(0);
1080   for(std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1081     {
1082       for(int j=0;j<(*fields)->getNumberOfFields();j++)
1083         {
1084           MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1085           std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1086           allFMTSLeavesToDisplaySafe.insert(allFMTSLeavesToDisplaySafe.end(),tmp.begin(),tmp.end());
1087         }
1088     }
1089   std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1090   for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1091     {
1092       allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1093     }
1094   std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1095   // memory safety part
1096   std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1097   for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1098     {
1099       allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1100       for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1101         {
1102           allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1103           allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1104         }
1105     }
1106   // end of memory safety part
1107   // 1st : timesteps, 2nd : meshName, 3rd : common support
1108   this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1109   for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1110     {
1111       vtksys_stl::vector< vtksys_stl::string > meshNamesLoc;
1112       std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1113       for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1114         {
1115           std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1116           vtksys_stl::vector< vtksys_stl::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1117           if(it==meshNamesLoc.end())
1118             {
1119               meshNamesLoc.push_back(meshName);
1120               splitByMeshName.resize(splitByMeshName.size()+1);
1121               splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1122             }
1123           else
1124             splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1125         }
1126       _data_structure[i].resize(meshNamesLoc.size());
1127       for(std::size_t j=0;j<splitByMeshName.size();j++)
1128         {
1129           std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFastCellSupportComparator> > fsp;
1130           std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1131           for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1132             sbmn[k]=splitByMeshName[j][k];
1133           //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1134           std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1135           std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1136           this->_data_structure[i][j].resize(commonSupSplit.size());
1137           for(std::size_t k=0;k<commonSupSplit.size();k++)
1138             {
1139               commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1140               for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1141                 {
1142                   commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1143                   commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1144                 }
1145             }
1146           for(std::size_t k=0;k<commonSupSplit.size();k++)
1147             this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1148         }
1149     }
1150   this->removeEmptyLeaves();
1151   this->assignIds();
1152 }
1153
1154 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1155 {
1156   std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1157   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1158     {
1159       std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1160       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1161         {
1162           std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1163           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1164             if(!(*it2).empty())
1165               newSD1.push_back(*it2);
1166           if(!newSD1.empty())
1167             newSD0.push_back(newSD1);
1168         }
1169       if(!newSD0.empty())
1170         newSD.push_back(newSD0);
1171     }
1172 }
1173
1174 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1175 {
1176   if(compInfos.size()!=1)
1177     return false;
1178   return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1179 }
1180
1181 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1182 {
1183   return _data_structure[0][0][0].getMeshName();
1184 }
1185
1186 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1187 {
1188   int lev1,lev2;
1189   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1190   return leaf.getTimeSteps(tk);
1191 }
1192
1193 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const
1194 {
1195   int lev0,lev1,lev2;
1196   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1197   meshName=leaf.getMeshName();
1198   std::vector<double> ts(leaf.getTimeSteps(tk));
1199   std::size_t zeTimeId(0);
1200   if(ts.size()!=1)
1201     {
1202       std::vector<double> ts2(ts.size());
1203       std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1204       std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1205       zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1206     }
1207   //2nd chance
1208   if(zeTimeId==(int)ts.size())
1209     zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1210   if(zeTimeId==(int)ts.size())
1211     {//OK the time requested does not fit time series given to ParaView. It is typically the case if more than one MEDReader instance are created or TimeInspector in real time mode.
1212       //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1213       int pos(-1);
1214       double valAttachedToPos(-std::numeric_limits<double>::max());
1215       for(std::size_t i=0;i<ts.size();i++)
1216         {
1217           if(ts[i]<timeReq)
1218             {
1219               if(ts[i]>valAttachedToPos)
1220                 {
1221                   pos=i;
1222                   valAttachedToPos=ts[i];
1223                 }
1224             }
1225         }
1226       if(pos==-1)
1227         {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1228           valAttachedToPos=std::numeric_limits<double>::max();
1229           for(std::size_t i=0;i<ts.size();i++)
1230             {
1231               if(ts[i]<valAttachedToPos)
1232                 {
1233                   pos=i;
1234                   valAttachedToPos=ts[i];
1235                 }
1236             }
1237         }
1238       zeTimeId=pos;
1239       std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1240       std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1241       oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1242       std::cerr << oss.str() << std::endl;
1243     }
1244   MEDTimeReq *tr(0);
1245   if(!isStdOrMode)
1246     tr=new MEDStdTimeReq((int)zeTimeId);
1247   else
1248     tr=new MEDModeTimeReq(tk.getTheVectOfBool());
1249   vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms));
1250   delete tr;
1251   return ret;
1252 }
1253
1254 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1255 {
1256   int nbOfActivated(0);
1257   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1258     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1259       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1260         if((*it2).isActivated())
1261           nbOfActivated++;
1262   if(nbOfActivated!=1)
1263     {
1264       std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1265       throw INTERP_KERNEL::Exception(oss.str().c_str());
1266     }
1267   int i0(0),i1(0),i2(0);
1268   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1269     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1270       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1271         if((*it2).isActivated())
1272           {
1273             lev0=i0; lev1=i1; lev2=i2;
1274             return *it2;
1275           }
1276   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1277 }
1278
1279 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1280 {
1281   int i(0);
1282   os << "#############################################" << std::endl;
1283   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1284     {
1285       int j(0);
1286       os << "TS" << i << std::endl;
1287       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1288         {
1289           int k(0);
1290           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1291             {
1292               if(k==0)
1293                 os << "   " << (*it2).getMeshName() << std::endl;
1294               os << "      Comp" << k  << std::endl;
1295               (*it2).printMySelf(os);
1296             }
1297         }
1298     }
1299     os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1300 }
1301
1302 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms, ParaMEDMEM::MEDFileFields *ret)
1303 {
1304   for(int i=0;i<ms->getNumberOfMeshes();i++)
1305     {
1306       MEDFileMesh *mm(ms->getMeshAtPos(i));
1307       MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1308       if(!mmu)
1309         continue;
1310       std::vector<int> levs(mm->getNonEmptyLevels());
1311       ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileField1TS> f1tsMultiLev(ParaMEDMEM::MEDFileField1TS::New());
1312       for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1313         {
1314           std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1315           for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1316             {
1317               ParaMEDMEM::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1318               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1319               f->setMesh(m);
1320               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1321               arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1322               arr->iota();
1323               f->setArray(arr);
1324               f->setName(mm->getName());
1325               f1tsMultiLev->setFieldNoProfileSBT(f);
1326             }
1327         }
1328       if(levs.empty())
1329         {
1330           std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1331           if(levsExt.size()==levs.size()+1)
1332             {
1333               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(1));
1334               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_NODES));
1335               f->setMesh(m);
1336               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1337               arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1338               arr->iota(); f->setArray(arr);
1339               f->setName(mm->getName());
1340               f1tsMultiLev->setFieldNoProfileSBT(f);
1341             }
1342           else
1343             continue;
1344         }
1345       //
1346       ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFieldMultiTS> fmtsMultiLev(ParaMEDMEM::MEDFileFieldMultiTS::New());
1347       fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1348       ret->pushField(fmtsMultiLev);
1349     }
1350 }
1351
1352 ParaMEDMEM::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms)
1353 {
1354   ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFields> ret(ParaMEDMEM::MEDFileFields::New());
1355   AppendFieldFromMeshes(ms,ret);
1356   return ret.retn();
1357 }
1358
1359 ///////////
1360
1361 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1362 {
1363 }
1364
1365 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1366 {
1367   switch(_policy)
1368     {
1369     case 0:
1370       return getTimeStepsRegardingPolicy0(tsPairs,ts);
1371     case 1:
1372       return getTimeStepsRegardingPolicy0(tsPairs,ts);
1373     default:
1374       throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1375     }
1376 }
1377
1378 /*!
1379  * policy = 0 :
1380  * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1381  * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1382  * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1383  */
1384 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1385 {
1386   std::size_t sz(ts.size());
1387   bool isInHumanRange(true);
1388   std::set<double> s;
1389   for(std::size_t i=0;i<sz;i++)
1390     {
1391       s.insert(ts[i]);
1392       if(ts[i]<=-1e299 || ts[i]>=1e299)
1393         isInHumanRange=false;
1394     }
1395   if(!isInHumanRange)
1396     return processedUsingPairOfIds(tsPairs);
1397   if(s.size()!=sz)
1398     return processedUsingPairOfIds(tsPairs);
1399   _postprocessed_time=ts;
1400   return getPostProcessedTime();
1401 }
1402
1403 /*!
1404  * policy = 1 :
1405  * idem than 0, except that ts is preaccumulated before invoking policy 0.
1406  */
1407 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1408 {
1409   std::size_t sz(ts.size());
1410   std::vector<double> ts2(sz);
1411   double acc(0.);
1412   for(std::size_t i=0;i<sz;i++)
1413     {
1414       ts2[i]=acc;
1415       acc+=ts[i];
1416     }
1417   return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1418 }
1419
1420 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1421 {
1422   std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1423   return (int)pos;
1424 }
1425
1426 void TimeKeeper::printSelf(std::ostream& oss) const
1427 {
1428   std::size_t sz(_activated_ts.size());
1429   for(std::size_t i=0;i<sz;i++)
1430     {
1431       oss << "(" << i << "," << _activated_ts[i].first << "), ";
1432     }
1433 }
1434
1435 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1436 {
1437   std::size_t sz(_activated_ts.size());
1438   std::vector<bool> ret(sz);
1439   for(std::size_t i=0;i<sz;i++)
1440     {
1441       ret[i]=_activated_ts[i].first;
1442     }
1443   return ret;
1444 }
1445
1446 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1447 {
1448   std::size_t sz(tsPairs.size());
1449   std::set<int> s0,s1;
1450   for(std::size_t i=0;i<sz;i++)
1451     { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1452   if(s0.size()==sz)
1453     {
1454       _postprocessed_time.resize(sz);
1455       for(std::size_t i=0;i<sz;i++)
1456         _postprocessed_time[i]=(double)tsPairs[i].first;
1457       return getPostProcessedTime();
1458     }
1459   if(s1.size()==sz)
1460     {
1461       _postprocessed_time.resize(sz);
1462       for(std::size_t i=0;i<sz;i++)
1463         _postprocessed_time[i]=(double)tsPairs[i].second;
1464       return getPostProcessedTime();
1465     }
1466   //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1467   _postprocessed_time.resize(sz);
1468   for(std::size_t i=0;i<sz;i++)
1469     _postprocessed_time[i]=(double)i;
1470   return getPostProcessedTime();
1471 }
1472
1473 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1474 {
1475   _activated_ts.resize(maxNumberOfTS);
1476   for(int i=0;i<maxNumberOfTS;i++)
1477     {
1478       std::ostringstream oss; oss << "000" << i;
1479       _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());
1480     }
1481 }