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