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