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