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