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