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