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