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