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