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