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