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