]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/MEDReader/IO/MEDFileFieldRepresentationTree.cxx
Salome HOME
new MEDReader
[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 void MEDFileFieldRepresentationLeaves::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
467 {
468   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
469     (*it).feedSIL(sil,root,edge,tsName,meshName,comSupStr,names);
470 }
471
472 bool MEDFileFieldRepresentationLeaves::containId(int id) const
473 {
474   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
475     if((*it).getId()==id)
476       return true;
477   return false;
478 }
479
480 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
481 {
482   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
483     if((*it).getZeName()==name)
484       {
485         id=(*it).getId();
486         return true;
487       }
488   return false;
489 }
490
491 bool MEDFileFieldRepresentationLeaves::isActivated() const
492 {
493   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
494     if((*it).getStatus())
495       return true;
496   return false;
497 }
498
499 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
500 {
501   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
502     {
503       os << "         - " << (*it0).getZeName() << " (";
504       if((*it0).getStatus())
505         os << "X";
506       else
507         os << " ";
508       os << ")" << std::endl;
509     }
510 }
511
512 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
513 {
514   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
515     (*it).setStatus(true);
516 }
517
518 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
519 {
520   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
521     if((*it).getId()==id)
522       return *it;
523   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
524 }
525
526 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
527 {
528   if(_arrays.size()<1)
529     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
530   std::vector<double> ret;
531   std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
532   return tk.getTimeStepsRegardingPolicy(dtits,ret);
533 }
534
535 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
536 {
537   if(!_arrays.empty())
538     return _arrays[0]->getTimeSteps(ts);
539   else
540     {
541       ts.clear();
542       return std::vector< std::pair<int,int> >();
543     }
544 }
545
546 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
547 {
548   std::ostringstream oss;
549   oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
550   std::vector<double> ret1;
551   std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
552   std::size_t sz(ret1.size());
553   for(std::size_t i=0;i<sz;i++)
554     {
555       oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
556       if(i!=sz-1)
557         oss << ", ";
558       std::string tmp(oss.str());
559       if(tmp.size()>200 && i!=sz-1)
560         {
561           oss << "...";
562           break;
563         }
564     }
565   oss << ")";
566   return oss.str();
567 }
568
569 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshes *meshes, vtkDataSet *ds) const
570 {
571   if(_arrays.size()<1)
572     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
573   MEDCouplingAutoRefCountObjectPtr<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
574   for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
575     if((*it).getStatus())
576       {
577         (*it).appendFields(tr,globs,mml,mst,ds);
578         (*it).appendELGAIfAny(ds);
579       }
580 }
581
582 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
583 {
584   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
585   DataArrayDouble *coordsMC(0);
586   DataArrayByte *typesMC(0);
587   DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
588   bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
589   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coordsSafe(coordsMC);
590   MEDCouplingAutoRefCountObjectPtr<DataArrayByte> typesSafe(typesMC);
591   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
592   //
593   int nbOfCells(typesSafe->getNbOfElems());
594   vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
595   vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
596   cellTypes->SetArray(reinterpret_cast<unsigned char *>(typesSafe->getPointer()),nbOfCells,0,VTK_DATA_ARRAY_FREE); typesSafe->accessToMemArray().setSpecificDeallocator(0);
597   vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
598   cellLocations->SetArray(cellLocationsSafe->getPointer(),nbOfCells,0,VTK_DATA_ARRAY_FREE); cellLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
599   vtkCellArray *cells(vtkCellArray::New());
600   vtkIdTypeArray *cells2(vtkIdTypeArray::New());
601   cells2->SetArray(cellsSafe->getPointer(),cellsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); cellsSafe->accessToMemArray().setSpecificDeallocator(0);
602   cells->SetCells(nbOfCells,cells2);
603   cells2->Delete();
604   if(faceLocationsMC!=0 && facesMC!=0)
605     {
606       vtkIdTypeArray *faces(vtkIdTypeArray::New());
607       faces->SetArray(facesSafe->getPointer(),facesSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); facesSafe->accessToMemArray().setSpecificDeallocator(0);
608       vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
609       faceLocations->SetArray(faceLocationsSafe->getPointer(),faceLocationsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); faceLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
610       ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
611       faceLocations->Delete();
612       faces->Delete();
613     }
614   else
615     ret->SetCells(cellTypes,cellLocations,cells);
616   cellTypes->Delete();
617   cellLocations->Delete();
618   cells->Delete();
619   vtkPoints *pts(vtkPoints::New());
620   vtkDoubleArray *pts2(vtkDoubleArray::New());
621   pts2->SetNumberOfComponents(3);
622   if(!statusOfCoords)
623     {
624       pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);
625       coordsSafe->accessToMemArray().setSpecificDeallocator(0);
626     }
627   else
628     pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
629   pts->SetData(pts2);
630   pts2->Delete();
631   ret->SetPoints(pts);
632   pts->Delete();
633   //
634   return ret;
635 }
636
637 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(ParaMEDMEM::MEDCMeshMultiLev *mm) const
638 {
639   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
640   bool isInternal;
641   std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
642   vtkDoubleArray *vtkTmp(0);
643   vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
644   std::size_t dim(arrs.size());
645   if(dim<1 || dim>3)
646     throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
647   int sizePerAxe[3]={1,1,1};
648   sizePerAxe[0]=arrs[0]->getNbOfElems();
649   if(dim>=2)
650     sizePerAxe[1]=arrs[1]->getNbOfElems();
651   if(dim==3)
652     sizePerAxe[2]=arrs[2]->getNbOfElems();
653   ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
654   vtkTmp=vtkDoubleArray::New();
655   vtkTmp->SetNumberOfComponents(1);
656   if(isInternal)
657     vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
658   else
659     { vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[0]->accessToMemArray().setSpecificDeallocator(0); }
660   ret->SetXCoordinates(vtkTmp);
661   vtkTmp->Delete();
662   arrs[0]->decrRef();
663   if(dim>=2)
664     {
665       vtkTmp=vtkDoubleArray::New();
666       vtkTmp->SetNumberOfComponents(1);
667       if(isInternal)
668         vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
669       else
670         { vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[1]->accessToMemArray().setSpecificDeallocator(0); }
671       ret->SetYCoordinates(vtkTmp);
672       vtkTmp->Delete();
673       arrs[1]->decrRef();
674     }
675   if(dim==3)
676     {
677       vtkTmp=vtkDoubleArray::New();
678       vtkTmp->SetNumberOfComponents(1);
679       if(isInternal)
680         vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
681       else
682         { vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[2]->accessToMemArray().setSpecificDeallocator(0); }
683       ret->SetZCoordinates(vtkTmp);
684       vtkTmp->Delete();
685       arrs[2]->decrRef();
686     }
687   return ret;
688 }
689
690 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(ParaMEDMEM::MEDCurveLinearMeshMultiLev *mm) const
691 {
692   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
693   int meshStr[3]={1,1,1};
694   DataArrayDouble *coords(0);
695   std::vector<int> nodeStrct;
696   bool isInternal;
697   mm->buildVTUArrays(coords,nodeStrct,isInternal);
698   std::size_t dim(nodeStrct.size());
699   if(dim<1 || dim>3)
700     throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
701   meshStr[0]=nodeStrct[0];
702   if(dim>=2)
703     meshStr[1]=nodeStrct[1];
704   if(dim==3)
705     meshStr[2]=nodeStrct[2];
706   vtkStructuredGrid *ret(vtkStructuredGrid::New());
707   ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
708   vtkDoubleArray *da(vtkDoubleArray::New());
709   da->SetNumberOfComponents(3);
710   if(coords->getNumberOfComponents()==3)
711     {
712       if(isInternal)
713         da->SetArray(coords->getPointer(),coords->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);//VTK has not the ownership of double * because MEDLoader main struct has it !
714       else
715         { da->SetArray(coords->getPointer(),coords->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); coords->accessToMemArray().setSpecificDeallocator(0); }
716     }
717   else
718     {
719       MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
720       da->SetArray(coords2->getPointer(),coords2->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);//let VTK deal with double *
721       coords2->accessToMemArray().setSpecificDeallocator(0);
722     }
723   coords->decrRef();
724   vtkPoints *points=vtkPoints::New();
725   ret->SetPoints(points);
726   points->SetData(da);
727   points->Delete();
728   da->Delete();
729   return ret;
730 }
731  
732 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDFileMeshes *meshes) const
733 {
734   static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
735   vtkDataSet *ret(0);
736   //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
737   MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change 
738   MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml2(mml->prepare());
739   MEDMeshMultiLev *ptMML2(mml2);
740   if(!_cached_ds)
741     {
742       MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
743       MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
744       MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
745       
746       if(ptUMML2)
747         {
748           ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
749         }
750       else if(ptCMML2)
751         {
752           ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
753         }
754       else if(ptCLMML2)
755         {
756           ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
757         }
758       else
759         throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
760       _cached_ds=ret->NewInstance();
761       _cached_ds->ShallowCopy(ret);
762     }
763   else
764     {
765       ret=_cached_ds->NewInstance();
766       ret->ShallowCopy(_cached_ds);
767     }
768   //
769   appendFields(tr,globs,mml,meshes,ret);
770   // The arrays links to mesh
771   DataArrayInt *famCells(0),*numCells(0);
772   bool noCpyFamCells(false),noCpyNumCells(false);
773   ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
774   if(famCells)
775     {
776       vtkIntArray *vtkTab(vtkIntArray::New());
777       vtkTab->SetNumberOfComponents(1);
778       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
779       if(noCpyFamCells)
780         vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
781       else
782         { vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famCells->accessToMemArray().setSpecificDeallocator(0); }
783       ret->GetCellData()->AddArray(vtkTab);
784       vtkTab->Delete();
785       famCells->decrRef();
786     }
787   ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
788   if(numCells)
789     {
790       vtkIntArray *vtkTab(vtkIntArray::New());
791       vtkTab->SetNumberOfComponents(1);
792       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
793       if(noCpyNumCells)
794         vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
795       else
796         { vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numCells->accessToMemArray().setSpecificDeallocator(0); }
797       ret->GetCellData()->AddArray(vtkTab);
798       vtkTab->Delete();
799       numCells->decrRef();
800     }
801   // The arrays links to mesh
802   DataArrayInt *famNodes(0),*numNodes(0);
803   bool noCpyFamNodes(false),noCpyNumNodes(false);
804   ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
805   if(famNodes)
806     {
807       vtkIntArray *vtkTab(vtkIntArray::New());
808       vtkTab->SetNumberOfComponents(1);
809       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
810       if(noCpyFamNodes)
811         vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
812       else
813         { vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famNodes->accessToMemArray().setSpecificDeallocator(0); }
814       ret->GetPointData()->AddArray(vtkTab);
815       vtkTab->Delete();
816       famNodes->decrRef();
817     }
818   ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
819   if(numNodes)
820     {
821       vtkIntArray *vtkTab(vtkIntArray::New());
822       vtkTab->SetNumberOfComponents(1);
823       vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
824       if(noCpyNumNodes)
825         vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
826       else
827         { vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numNodes->accessToMemArray().setSpecificDeallocator(0); }
828       ret->GetPointData()->AddArray(vtkTab);
829       vtkTab->Delete();
830       numNodes->decrRef();
831     }
832   return ret;
833 }
834
835 //////////////////////
836
837 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
838 {
839 }
840
841 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
842 {
843   int ret(0);
844   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
845     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
846       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
847         ret+=(*it2).getNumberOfArrays();
848   return ret;
849 }
850
851 void MEDFileFieldRepresentationTree::assignIds() const
852 {
853   int zeId(0);
854   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
855     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
856       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
857         (*it2).setId(zeId);
858 }
859 void MEDFileFieldRepresentationTree::activateTheFirst() const
860 {
861   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
862     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
863       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
864         {
865           (*it2).activateAllArrays();
866           return ;
867         }
868 }
869
870 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
871 {
872   std::size_t it0Cnt(0);
873   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
874     {
875       vtkIdType InfoOnTSId(sil->AddChild(root,edge));
876       names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
877       //
878       vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
879       std::vector<double> ts;
880       std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
881       std::size_t nbOfTS(dtits.size());
882       std::ostringstream oss3; oss3 << nbOfTS;
883       names.push_back(oss3.str());
884       for(std::size_t i=0;i<nbOfTS;i++)
885         {
886           std::ostringstream oss4; oss4 << dtits[i].first;
887           vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
888           names.push_back(oss4.str());
889           std::ostringstream oss5; oss5 << dtits[i].second;
890           vtkIdType ItId(sil->AddChild(DtId,edge));
891           names.push_back(oss5.str());
892           std::ostringstream oss6; oss6 << ts[i];
893           sil->AddChild(ItId,edge);
894           names.push_back(oss6.str());
895         }
896       //
897       std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
898       std::string tsName(oss.str());
899       vtkIdType typeId0(sil->AddChild(root,edge));
900       names.push_back(tsName);
901       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
902         {
903           std::string meshName((*it1)[0].getMeshName());
904           vtkIdType typeId1(sil->AddChild(typeId0,edge));
905           names.push_back(meshName);
906           std::size_t it2Cnt(0);
907           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
908             {
909               std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
910               std::string comSupStr(oss2.str());
911               vtkIdType typeId2(sil->AddChild(typeId1,edge));
912               names.push_back(comSupStr);
913               (*it2).feedSIL(sil,typeId2,edge,tsName,meshName,comSupStr,names);
914             } 
915         }
916     }
917 }
918
919 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
920 {
921   int dummy0(0),dummy1(0),dummy2(0);
922   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
923   std::string ret(leaf.getMeshName());
924   int i(0);
925   MEDFileMesh *m(0);
926   for(;i<_ms->getNumberOfMeshes();i++)
927     {
928       m=_ms->getMeshAtPos(i);
929       if(m->getName()==ret)
930         break;
931     }
932   if(i==_ms->getNumberOfMeshes())
933     throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
934   vtkIdType typeId0(sil->AddChild(root,edge));
935   names.push_back(m->getName());
936   //
937   vtkIdType typeId1(sil->AddChild(typeId0,edge));
938   names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
939   std::vector<std::string> grps(m->getGroupsNames());
940   for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
941     {
942       vtkIdType typeId2(sil->AddChild(typeId1,edge));
943       names.push_back(*it0);
944       std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
945       for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
946         {
947           sil->AddChild(typeId2,edge);
948           names.push_back((*it1).c_str());
949         }
950     }
951   //
952   vtkIdType typeId11(sil->AddChild(typeId0,edge));
953   names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
954   std::vector<std::string> fams(m->getFamiliesNames());
955   for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
956     {
957       sil->AddChild(typeId11,edge);
958       int famId(m->getFamilyId((*it00).c_str()));
959       std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
960       names.push_back(oss.str());
961     }
962   return ret;
963 }
964
965 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
966 {
967   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
968     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
969       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
970         if((*it2).containId(id))
971           return (*it2).getLeafArr(id);
972   throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
973 }
974
975 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
976 {
977   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
978   return elt.getZeName();
979 }
980
981 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
982 {
983   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
984   return elt.getStatus();
985 }
986
987 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
988 {
989   int ret(-1);
990   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
991     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
992       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
993         if((*it2).containZeName(name,ret))
994           return ret;
995   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name !");
996 }
997
998 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
999 {
1000   const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1001   bool ret(elt.setStatus(status));//to be implemented
1002   return ret;
1003 }
1004
1005 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1006 {
1007   int ret(0);
1008   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1009     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1010       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1011         ret=std::max(ret,(*it2).getNumberOfTS());
1012   return ret;
1013 }
1014
1015 /*!
1016  * 
1017  */
1018 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv)
1019 {
1020   if(isMEDOrSauv)
1021     {
1022       _ms=MEDFileMeshes::New(fileName);
1023       _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
1024     }
1025   else
1026     {
1027       MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::SauvReader> sr(ParaMEDMEM::SauvReader::New(fileName));
1028       MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileData> mfd(sr->loadInMEDFileDS());
1029       _ms=mfd->getMeshes(); _ms->incrRef();
1030       int nbMeshes(_ms->getNumberOfMeshes());
1031       for(int i=0;i<nbMeshes;i++)
1032         {
1033           ParaMEDMEM::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1034           ParaMEDMEM::MEDFileUMesh *tmp2(dynamic_cast<ParaMEDMEM::MEDFileUMesh *>(tmp));
1035           if(tmp2)
1036             tmp2->forceComputationOfParts();
1037         }
1038       _fields=mfd->getFields();
1039       if((ParaMEDMEM::MEDFileFields *)_fields)
1040         _fields->incrRef();
1041     }
1042   if(!((ParaMEDMEM::MEDFileFields *)_fields))
1043     {
1044       _fields=BuildFieldFromMeshes(_ms);
1045     }
1046   else
1047     {
1048       AppendFieldFromMeshes(_ms,_fields);
1049     }
1050   std::vector<std::string> meshNames(_ms->getMeshesNames());
1051   std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> > fields_per_mesh(meshNames.size());
1052   for(std::size_t i=0;i<meshNames.size();i++)
1053     {
1054       fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1055     }
1056   std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1057   std::size_t k(0);
1058   for(std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1059     {
1060       for(int j=0;j<(*fields)->getNumberOfFields();j++)
1061         {
1062           MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1063           std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1064           allFMTSLeavesToDisplaySafe.insert(allFMTSLeavesToDisplaySafe.end(),tmp.begin(),tmp.end());
1065         }
1066     }
1067   std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1068   for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1069     {
1070       allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1071     }
1072   std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1073   // memory safety part
1074   std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1075   for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1076     {
1077       allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1078       for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1079         {
1080           allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1081           allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1082         }
1083     }
1084   // end of memory safety part
1085   // 1st : timesteps, 2nd : meshName, 3rd : common support
1086   this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1087   for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1088     {
1089       vtksys_stl::vector< vtksys_stl::string > meshNamesLoc;
1090       std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1091       for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1092         {
1093           std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1094           vtksys_stl::vector< vtksys_stl::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1095           if(it==meshNamesLoc.end())
1096             {
1097               meshNamesLoc.push_back(meshName);
1098               splitByMeshName.resize(splitByMeshName.size()+1);
1099               splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1100             }
1101           else
1102             splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1103         }
1104       _data_structure[i].resize(meshNamesLoc.size());
1105       for(std::size_t j=0;j<splitByMeshName.size();j++)
1106         {
1107           std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFastCellSupportComparator> > fsp;
1108           std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1109           for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1110             sbmn[k]=splitByMeshName[j][k];
1111           //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1112           std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1113           std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1114           this->_data_structure[i][j].resize(commonSupSplit.size());
1115           for(std::size_t k=0;k<commonSupSplit.size();k++)
1116             {
1117               commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1118               for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1119                 {
1120                   commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1121                   commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1122                 }
1123             }
1124           for(std::size_t k=0;k<commonSupSplit.size();k++)
1125             this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1126         }
1127     }
1128   this->removeEmptyLeaves();
1129   this->assignIds();
1130 }
1131
1132 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1133 {
1134   std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1135   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1136     {
1137       std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1138       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1139         {
1140           std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1141           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1142             if(!(*it2).empty())
1143               newSD1.push_back(*it2);
1144           if(!newSD1.empty())
1145             newSD0.push_back(newSD1);
1146         }
1147       if(!newSD0.empty())
1148         newSD.push_back(newSD0);
1149     }
1150 }
1151
1152 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1153 {
1154   if(compInfos.size()!=1)
1155     return false;
1156   return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1157 }
1158
1159 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1160 {
1161   return _data_structure[0][0][0].getMeshName();
1162 }
1163
1164 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1165 {
1166   int lev1,lev2;
1167   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1168   return leaf.getTimeSteps(tk);
1169 }
1170
1171 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const
1172 {
1173   int lev0,lev1,lev2;
1174   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1175   meshName=leaf.getMeshName();
1176   std::vector<double> ts(leaf.getTimeSteps(tk));
1177   std::size_t zeTimeId(0);
1178   if(ts.size()!=1)
1179     {
1180       std::vector<double> ts2(ts.size());
1181       std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1182       std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1183       zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1184     }
1185   //2nd chance
1186   if(zeTimeId==(int)ts.size())
1187     zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1188   if(zeTimeId==(int)ts.size())
1189     {//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.
1190       //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1191       int pos(-1);
1192       double valAttachedToPos(-std::numeric_limits<double>::max());
1193       for(std::size_t i=0;i<ts.size();i++)
1194         {
1195           if(ts[i]<timeReq)
1196             {
1197               if(ts[i]>valAttachedToPos)
1198                 {
1199                   pos=i;
1200                   valAttachedToPos=ts[i];
1201                 }
1202             }
1203         }
1204       if(pos==-1)
1205         {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1206           valAttachedToPos=std::numeric_limits<double>::max();
1207           for(std::size_t i=0;i<ts.size();i++)
1208             {
1209               if(ts[i]<valAttachedToPos)
1210                 {
1211                   pos=i;
1212                   valAttachedToPos=ts[i];
1213                 }
1214             }
1215         }
1216       zeTimeId=pos;
1217       std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1218       std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1219       oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1220       std::cerr << oss.str() << std::endl;
1221     }
1222   MEDTimeReq *tr(0);
1223   if(!isStdOrMode)
1224     tr=new MEDStdTimeReq((int)zeTimeId);
1225   else
1226     tr=new MEDModeTimeReq(tk.getTheVectOfBool());
1227   vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms));
1228   delete tr;
1229   return ret;
1230 }
1231
1232 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1233 {
1234   int nbOfActivated(0);
1235   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1236     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1237       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1238         if((*it2).isActivated())
1239           nbOfActivated++;
1240   if(nbOfActivated!=1)
1241     {
1242       std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1243       throw INTERP_KERNEL::Exception(oss.str().c_str());
1244     }
1245   int i0(0),i1(0),i2(0);
1246   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1247     for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1248       for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1249         if((*it2).isActivated())
1250           {
1251             lev0=i0; lev1=i1; lev2=i2;
1252             return *it2;
1253           }
1254   throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1255 }
1256
1257 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1258 {
1259   int i(0);
1260   os << "#############################################" << std::endl;
1261   for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1262     {
1263       int j(0);
1264       os << "TS" << i << std::endl;
1265       for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1266         {
1267           int k(0);
1268           for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1269             {
1270               if(k==0)
1271                 os << "   " << (*it2).getMeshName() << std::endl;
1272               os << "      Comp" << k  << std::endl;
1273               (*it2).printMySelf(os);
1274             }
1275         }
1276     }
1277     os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1278 }
1279
1280 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms, ParaMEDMEM::MEDFileFields *ret)
1281 {
1282   for(int i=0;i<ms->getNumberOfMeshes();i++)
1283     {
1284       MEDFileMesh *mm(ms->getMeshAtPos(i));
1285       MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1286       if(!mmu)
1287         continue;
1288       std::vector<int> levs(mm->getNonEmptyLevels());
1289       ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileField1TS> f1tsMultiLev(ParaMEDMEM::MEDFileField1TS::New());
1290       for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1291         {
1292           std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1293           for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1294             {
1295               ParaMEDMEM::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1296               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1297               f->setMesh(m);
1298               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1299               arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1300               arr->iota();
1301               f->setArray(arr);
1302               f->setName(mm->getName());
1303               f1tsMultiLev->setFieldNoProfileSBT(f);
1304             }
1305         }
1306       if(levs.empty())
1307         {
1308           std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1309           if(levsExt.size()==levs.size()+1)
1310             {
1311               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(1));
1312               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_NODES));
1313               f->setMesh(m);
1314               ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1315               arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1316               arr->iota(); f->setArray(arr);
1317               f->setName(mm->getName());
1318               f1tsMultiLev->setFieldNoProfileSBT(f);
1319             }
1320           else
1321             continue;
1322         }
1323       //
1324       ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFieldMultiTS> fmtsMultiLev(ParaMEDMEM::MEDFileFieldMultiTS::New());
1325       fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1326       ret->pushField(fmtsMultiLev);
1327     }
1328 }
1329
1330 ParaMEDMEM::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms)
1331 {
1332   ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFields> ret(ParaMEDMEM::MEDFileFields::New());
1333   AppendFieldFromMeshes(ms,ret);
1334   return ret.retn();
1335 }
1336
1337 ///////////
1338
1339 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1340 {
1341 }
1342
1343 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1344 {
1345   switch(_policy)
1346     {
1347     case 0:
1348       return getTimeStepsRegardingPolicy0(tsPairs,ts);
1349     case 1:
1350       return getTimeStepsRegardingPolicy0(tsPairs,ts);
1351     default:
1352       throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1353     }
1354 }
1355
1356 /*!
1357  * policy = 0 :
1358  * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1359  * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1360  * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1361  */
1362 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1363 {
1364   std::size_t sz(ts.size());
1365   bool isInHumanRange(true);
1366   std::set<double> s;
1367   for(std::size_t i=0;i<sz;i++)
1368     {
1369       s.insert(ts[i]);
1370       if(ts[i]<=-1e299 || ts[i]>=1e299)
1371         isInHumanRange=false;
1372     }
1373   if(!isInHumanRange)
1374     return processedUsingPairOfIds(tsPairs);
1375   if(s.size()!=sz)
1376     return processedUsingPairOfIds(tsPairs);
1377   _postprocessed_time=ts;
1378   return getPostProcessedTime();
1379 }
1380
1381 /*!
1382  * policy = 1 :
1383  * idem than 0, except that ts is preaccumulated before invoking policy 0.
1384  */
1385 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1386 {
1387   std::size_t sz(ts.size());
1388   std::vector<double> ts2(sz);
1389   double acc(0.);
1390   for(std::size_t i=0;i<sz;i++)
1391     {
1392       ts2[i]=acc;
1393       acc+=ts[i];
1394     }
1395   return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1396 }
1397
1398 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1399 {
1400   std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1401   return (int)pos;
1402 }
1403
1404 void TimeKeeper::printSelf(std::ostream& oss) const
1405 {
1406   std::size_t sz(_activated_ts.size());
1407   for(std::size_t i=0;i<sz;i++)
1408     {
1409       oss << "(" << i << "," << _activated_ts[i].first << "), ";
1410     }
1411 }
1412
1413 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1414 {
1415   std::size_t sz(_activated_ts.size());
1416   std::vector<bool> ret(sz);
1417   for(std::size_t i=0;i<sz;i++)
1418     {
1419       ret[i]=_activated_ts[i].first;
1420     }
1421   return ret;
1422 }
1423
1424 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1425 {
1426   std::size_t sz(tsPairs.size());
1427   std::set<int> s0,s1;
1428   for(std::size_t i=0;i<sz;i++)
1429     { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1430   if(s0.size()==sz)
1431     {
1432       _postprocessed_time.resize(sz);
1433       for(std::size_t i=0;i<sz;i++)
1434         _postprocessed_time[i]=(double)tsPairs[i].first;
1435       return getPostProcessedTime();
1436     }
1437   if(s1.size()==sz)
1438     {
1439       _postprocessed_time.resize(sz);
1440       for(std::size_t i=0;i<sz;i++)
1441         _postprocessed_time[i]=(double)tsPairs[i].second;
1442       return getPostProcessedTime();
1443     }
1444   //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1445   _postprocessed_time.resize(sz);
1446   for(std::size_t i=0;i<sz;i++)
1447     _postprocessed_time[i]=(double)i;
1448   return getPostProcessedTime();
1449 }
1450
1451 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1452 {
1453   _activated_ts.resize(maxNumberOfTS);
1454   for(int i=0;i<maxNumberOfTS;i++)
1455     {
1456       std::ostringstream oss; oss << "000" << i;
1457       _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());
1458     }
1459 }