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