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