1 // Copyright (C) 2010-2014 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay
21 #include "MEDTimeReq.hxx"
22 #include "MEDUtilities.hxx"
24 #include "MEDFileFieldRepresentationTree.hxx"
25 #include "MEDCouplingFieldDiscretization.hxx"
26 #include "MEDCouplingFieldDouble.hxx"
27 #include "InterpKernelGaussCoords.hxx"
28 #include "MEDFileData.hxx"
29 #include "SauvReader.hxx"
31 #include "vtkXMLUnstructuredGridWriter.h"//
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"
49 #include "vtksys/stl/string"
50 #include "vtksys/ios/fstream"
51 #include "vtksys/stl/algorithm"
52 #include "vtkMutableDirectedGraph.h"
54 using namespace ParaMEDMEM;
56 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
58 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
60 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
62 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
64 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
66 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
68 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
70 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
72 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
74 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
76 vtkIdTypeArray *ELGACmp::findOrCreate(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, bool& isNew) const
78 vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
87 return createNew(globs,locsReallyUsed,vtkd,ds);
91 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd) const
93 std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
94 if(it==_loc_names.end())
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++)
101 key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
103 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
107 vtkIdTypeArray *ELGACmp::createNew(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds) const
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;
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++)
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++)
138 const double *pt0(calculator.getFunctionValues(i));
139 if(ct!=INTERP_KERNEL::NORM_HEXA27)
140 std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
143 for(int j=0;j<27;j++)
144 shape[nbPtsPerCell*i+j]=pt0[MEDMeshMultiLev::HEXA27_PERM_ARRAY[j]];
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]));
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));
156 vtkIdType ncell(ds->GetNumberOfCells());
157 int *pt(new int[ncell]),offset(0);
158 for(vtkIdType cellId=0;cellId<ncell;cellId++)
160 vtkCell *cell(ds->GetCell(cellId));
161 int delta(m[cell->GetCellType()]);
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);
176 _defs.push_back(defs);
179 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
181 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
182 ds->GetCellData()->AddArray(*it);
187 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
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();
196 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
200 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>& arr):ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
202 std::vector< std::vector<ParaMEDMEM::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
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();
212 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
214 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>::operator=(other);
217 _ze_name=other._ze_name;
218 _ze_full_name.clear();
222 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
227 int MEDFileFieldRepresentationLeavesArrays::getId() const
232 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
234 return _ze_full_name;
237 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
239 vtkIdType refId(sil->AddChild(root,edge));
240 names.push_back(_ze_name);
242 if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
244 sil->AddChild(refId,edge);
245 names.push_back(std::string());
249 void MEDFileFieldRepresentationLeavesArrays::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
251 std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
252 _ze_full_name=oss3.str();
255 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
260 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
262 bool ret(_activated!=status);
267 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshStruct *mst, vtkDataSet *ds) const
269 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
270 static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
271 tr->setNumberOfTS((operator->())->getNumberOfTS());
273 for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
275 MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeField1TS> f1ts((operator->())->getTimeStepAtPos(tr->getCurrent()));
276 MEDFileAnyTypeField1TS *f1tsPtr(f1ts);
277 MEDFileField1TS *f1tsPtrDbl(dynamic_cast<MEDFileField1TS *>(f1tsPtr));
278 MEDFileIntField1TS *f1tsPtrInt(dynamic_cast<MEDFileIntField1TS *>(f1tsPtr));
279 DataArray *crudeArr(0),*postProcessedArr(0);
281 crudeArr=f1tsPtrDbl->getUndergroundDataArray();
283 crudeArr=f1tsPtrInt->getUndergroundDataArray();
285 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment !");
286 MEDFileField1TSStructItem fsst(MEDFileField1TSStructItem::BuildItemFrom(f1ts,mst));
287 f1ts->loadArraysIfNecessary();
288 MEDCouplingAutoRefCountObjectPtr<DataArray> v(mml->buildDataArray(fsst,globs,crudeArr));
291 std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
293 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
294 vtkFieldData *att(0);
299 att=ds->GetCellData();
304 att=ds->GetPointData();
309 att=ds->GetFieldData();
314 att=ds->GetFieldData();
318 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
322 DataArray *vPtr(v); DataArrayDouble *vd(static_cast<DataArrayDouble *>(vPtr));
323 vtkDoubleArray *vtkd(vtkDoubleArray::New());
324 vtkd->SetNumberOfComponents(vd->getNumberOfComponents());
325 for(int i=0;i<vd->getNumberOfComponents();i++)
326 vtkd->SetComponentName(i,vd->getInfoOnComponent(i).c_str());
327 if(postProcessedArr!=crudeArr)
329 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vd->accessToMemArray().setSpecificDeallocator(0);
333 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
335 std::string name(tr->buildName(f1ts->getName()));
336 vtkd->SetName(name.c_str());
339 if(discs[0]==ON_GAUSS_PT)
342 _elga_cmp.findOrCreate(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp);
344 if(discs[0]==ON_GAUSS_NE)
346 vtkIdTypeArray *elno(vtkIdTypeArray::New());
347 elno->SetNumberOfComponents(1);
348 vtkIdType ncell(ds->GetNumberOfCells());
349 int *pt(new int[ncell]),offset(0);
350 std::set<int> cellTypes;
351 for(vtkIdType cellId=0;cellId<ncell;cellId++)
353 vtkCell *cell(ds->GetCell(cellId));
354 int delta(cell->GetNumberOfPoints());
355 cellTypes.insert(cell->GetCellType());
359 elno->GetInformation()->Set(MEDUtilities::ELNO(),1);
360 elno->SetArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
361 std::string nameElno("ELNO"); nameElno+="@"; nameElno+=name;
362 elno->SetName(nameElno.c_str());
363 ds->GetCellData()->AddArray(elno);
364 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elno->GetName());
365 elno->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
367 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
368 for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
370 const unsigned char *pos(std::find(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,*it));
371 if(pos==MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
373 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)std::distance(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos));
374 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
375 int nbGaussPt(cm.getNumberOfNodes()),dim(cm.getDimension());
376 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
377 double *shape(new double[nbGaussPt*nbGaussPt]);
379 const double *gsCoords(MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType(ct,dummy));
380 const double *refCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));
381 const double *weights(MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(ct,dummy));
382 std::vector<double> gsCoords2(gsCoords,gsCoords+nbGaussPt*dim),refCoords2(refCoords,refCoords+nbGaussPt*dim);
383 INTERP_KERNEL::GaussInfo calculator(ct,gsCoords2,nbGaussPt,refCoords2,nbGaussPt);
384 calculator.initLocalInfo();
385 for(int i=0;i<nbGaussPt;i++)
387 const double *pt0(calculator.getFunctionValues(i));
388 std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
390 def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
392 key->Set(elno->GetInformation(),def,*it);
393 key->Set(vtkd->GetInformation(),def,*it);
402 DataArray *vPtr(v); DataArrayInt *vi(static_cast<DataArrayInt *>(vPtr));
403 vtkIntArray *vtkd(vtkIntArray::New());
404 vtkd->SetNumberOfComponents(vi->getNumberOfComponents());
405 for(int i=0;i<vi->getNumberOfComponents();i++)
406 vtkd->SetComponentName(i,vi->getVarOnComponent(i).c_str());
407 if(postProcessedArr!=crudeArr)
409 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vi->accessToMemArray().setSpecificDeallocator(0);
413 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
415 std::string name(tr->buildName(f1ts->getName()));
416 vtkd->SetName(name.c_str());
421 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
425 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
427 _elga_cmp.appendELGAIfAny(ds);
432 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
436 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves(const std::vector< ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS> >& arr,
437 const ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFastCellSupportComparator>& fsp):_arrays(arr.size()),_fsp(fsp),_cached_ds(0)
439 for(std::size_t i=0;i<arr.size();i++)
440 _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
443 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
446 _cached_ds->Delete();
449 bool MEDFileFieldRepresentationLeaves::empty() const
451 const MEDFileFastCellSupportComparator *fcscp(_fsp);
452 return fcscp==0 || _arrays.empty();
455 void MEDFileFieldRepresentationLeaves::setId(int& id) const
457 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
461 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
463 return _arrays[0]->getMeshName();
466 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
468 return (int)_arrays.size();
471 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
473 return _arrays[0]->getNumberOfTS();
476 void MEDFileFieldRepresentationLeaves::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
478 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
479 (*it).computeFullNameInLeaves(tsName,meshName,comSupStr);
483 * \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.
485 void MEDFileFieldRepresentationLeaves::feedSIL(const ParaMEDMEM::MEDFileMeshes *ms, const std::string& meshName, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
487 vtkIdType root2(sil->AddChild(root,edge));
488 names.push_back(std::string("Arrs"));
489 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
490 (*it).feedSIL(sil,root2,edge,names);
492 vtkIdType root3(sil->AddChild(root,edge));
493 names.push_back(std::string("InfoOnGeoType"));
494 const ParaMEDMEM::MEDFileMesh *m(0);
496 m=ms->getMeshWithName(meshName);
497 const ParaMEDMEM::MEDFileFastCellSupportComparator *fsp(_fsp);
498 if(!fsp || fsp->getNumberOfTS()==0)
500 std::vector< INTERP_KERNEL::NormalizedCellType > gts(fsp->getGeoTypesAt(0,m));
501 for(std::vector< INTERP_KERNEL::NormalizedCellType >::const_iterator it2=gts.begin();it2!=gts.end();it2++)
503 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(*it2));
504 std::string cmStr(cm.getRepr()); cmStr=cmStr.substr(5);//skip "NORM_"
505 sil->AddChild(root3,edge);
506 names.push_back(cmStr);
510 bool MEDFileFieldRepresentationLeaves::containId(int id) const
512 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
513 if((*it).getId()==id)
518 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
520 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
521 if((*it).getZeName()==name)
529 bool MEDFileFieldRepresentationLeaves::isActivated() const
531 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
532 if((*it).getStatus())
537 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
539 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
541 os << " - " << (*it0).getZeName() << " (";
542 if((*it0).getStatus())
546 os << ")" << std::endl;
550 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
552 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
553 (*it).setStatus(true);
556 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
558 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
559 if((*it).getId()==id)
561 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
564 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
567 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
568 std::vector<double> ret;
569 std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
570 return tk.getTimeStepsRegardingPolicy(dtits,ret);
573 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
576 return _arrays[0]->getTimeSteps(ts);
580 return std::vector< std::pair<int,int> >();
584 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
586 std::ostringstream oss;
587 oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
588 std::vector<double> ret1;
589 std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
590 std::size_t sz(ret1.size());
591 for(std::size_t i=0;i<sz;i++)
593 oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
596 std::string tmp(oss.str());
597 if(tmp.size()>200 && i!=sz-1)
607 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshes *meshes, vtkDataSet *ds) const
610 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
611 MEDCouplingAutoRefCountObjectPtr<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
612 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
613 if((*it).getStatus())
615 (*it).appendFields(tr,globs,mml,mst,ds);
616 (*it).appendELGAIfAny(ds);
620 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
622 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
623 DataArrayDouble *coordsMC(0);
624 DataArrayByte *typesMC(0);
625 DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
626 bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
627 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coordsSafe(coordsMC);
628 MEDCouplingAutoRefCountObjectPtr<DataArrayByte> typesSafe(typesMC);
629 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
631 int nbOfCells(typesSafe->getNbOfElems());
632 vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
633 vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
634 cellTypes->SetArray(reinterpret_cast<unsigned char *>(typesSafe->getPointer()),nbOfCells,0,VTK_DATA_ARRAY_FREE); typesSafe->accessToMemArray().setSpecificDeallocator(0);
635 vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
636 cellLocations->SetArray(cellLocationsSafe->getPointer(),nbOfCells,0,VTK_DATA_ARRAY_FREE); cellLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
637 vtkCellArray *cells(vtkCellArray::New());
638 vtkIdTypeArray *cells2(vtkIdTypeArray::New());
639 cells2->SetArray(cellsSafe->getPointer(),cellsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); cellsSafe->accessToMemArray().setSpecificDeallocator(0);
640 cells->SetCells(nbOfCells,cells2);
642 if(faceLocationsMC!=0 && facesMC!=0)
644 vtkIdTypeArray *faces(vtkIdTypeArray::New());
645 faces->SetArray(facesSafe->getPointer(),facesSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); facesSafe->accessToMemArray().setSpecificDeallocator(0);
646 vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
647 faceLocations->SetArray(faceLocationsSafe->getPointer(),faceLocationsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); faceLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
648 ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
649 faceLocations->Delete();
653 ret->SetCells(cellTypes,cellLocations,cells);
655 cellLocations->Delete();
657 vtkPoints *pts(vtkPoints::New());
658 vtkDoubleArray *pts2(vtkDoubleArray::New());
659 pts2->SetNumberOfComponents(3);
662 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);
663 coordsSafe->accessToMemArray().setSpecificDeallocator(0);
666 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
675 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(ParaMEDMEM::MEDCMeshMultiLev *mm) const
677 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
679 std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
680 vtkDoubleArray *vtkTmp(0);
681 vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
682 std::size_t dim(arrs.size());
684 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
685 int sizePerAxe[3]={1,1,1};
686 sizePerAxe[0]=arrs[0]->getNbOfElems();
688 sizePerAxe[1]=arrs[1]->getNbOfElems();
690 sizePerAxe[2]=arrs[2]->getNbOfElems();
691 ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
692 vtkTmp=vtkDoubleArray::New();
693 vtkTmp->SetNumberOfComponents(1);
695 vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
697 { vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[0]->accessToMemArray().setSpecificDeallocator(0); }
698 ret->SetXCoordinates(vtkTmp);
703 vtkTmp=vtkDoubleArray::New();
704 vtkTmp->SetNumberOfComponents(1);
706 vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
708 { vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[1]->accessToMemArray().setSpecificDeallocator(0); }
709 ret->SetYCoordinates(vtkTmp);
715 vtkTmp=vtkDoubleArray::New();
716 vtkTmp->SetNumberOfComponents(1);
718 vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
720 { vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[2]->accessToMemArray().setSpecificDeallocator(0); }
721 ret->SetZCoordinates(vtkTmp);
728 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(ParaMEDMEM::MEDCurveLinearMeshMultiLev *mm) const
730 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
731 int meshStr[3]={1,1,1};
732 DataArrayDouble *coords(0);
733 std::vector<int> nodeStrct;
735 mm->buildVTUArrays(coords,nodeStrct,isInternal);
736 std::size_t dim(nodeStrct.size());
738 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
739 meshStr[0]=nodeStrct[0];
741 meshStr[1]=nodeStrct[1];
743 meshStr[2]=nodeStrct[2];
744 vtkStructuredGrid *ret(vtkStructuredGrid::New());
745 ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
746 vtkDoubleArray *da(vtkDoubleArray::New());
747 da->SetNumberOfComponents(3);
748 if(coords->getNumberOfComponents()==3)
751 da->SetArray(coords->getPointer(),coords->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);//VTK has not the ownership of double * because MEDLoader main struct has it !
753 { da->SetArray(coords->getPointer(),coords->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); coords->accessToMemArray().setSpecificDeallocator(0); }
757 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
758 da->SetArray(coords2->getPointer(),coords2->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);//let VTK deal with double *
759 coords2->accessToMemArray().setSpecificDeallocator(0);
762 vtkPoints *points=vtkPoints::New();
763 ret->SetPoints(points);
770 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDFileMeshes *meshes) const
772 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
774 //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
775 MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change
776 MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml2(mml->prepare());
777 MEDMeshMultiLev *ptMML2(mml2);
780 MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
781 MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
782 MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
786 ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
790 ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
794 ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
797 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
798 _cached_ds=ret->NewInstance();
799 _cached_ds->ShallowCopy(ret);
803 ret=_cached_ds->NewInstance();
804 ret->ShallowCopy(_cached_ds);
807 appendFields(tr,globs,mml,meshes,ret);
808 // The arrays links to mesh
809 DataArrayInt *famCells(0),*numCells(0);
810 bool noCpyFamCells(false),noCpyNumCells(false);
811 ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
814 vtkIntArray *vtkTab(vtkIntArray::New());
815 vtkTab->SetNumberOfComponents(1);
816 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
818 vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
820 { vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famCells->accessToMemArray().setSpecificDeallocator(0); }
821 ret->GetCellData()->AddArray(vtkTab);
825 ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
828 vtkIntArray *vtkTab(vtkIntArray::New());
829 vtkTab->SetNumberOfComponents(1);
830 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
832 vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
834 { vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numCells->accessToMemArray().setSpecificDeallocator(0); }
835 ret->GetCellData()->AddArray(vtkTab);
839 // The arrays links to mesh
840 DataArrayInt *famNodes(0),*numNodes(0);
841 bool noCpyFamNodes(false),noCpyNumNodes(false);
842 ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
845 vtkIntArray *vtkTab(vtkIntArray::New());
846 vtkTab->SetNumberOfComponents(1);
847 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
849 vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
851 { vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famNodes->accessToMemArray().setSpecificDeallocator(0); }
852 ret->GetPointData()->AddArray(vtkTab);
856 ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
859 vtkIntArray *vtkTab(vtkIntArray::New());
860 vtkTab->SetNumberOfComponents(1);
861 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
863 vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
865 { vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numNodes->accessToMemArray().setSpecificDeallocator(0); }
866 ret->GetPointData()->AddArray(vtkTab);
873 //////////////////////
875 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
879 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
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 ret+=(*it2).getNumberOfArrays();
889 void MEDFileFieldRepresentationTree::assignIds() const
892 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
893 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
894 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
898 void MEDFileFieldRepresentationTree::computeFullNameInLeaves() const
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++)
903 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
904 std::string tsName(oss.str());
905 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
907 std::string meshName((*it1)[0].getMeshName());
908 std::size_t it2Cnt(0);
909 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
911 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
912 std::string comSupStr(oss2.str());
913 (*it2).computeFullNameInLeaves(tsName,meshName,comSupStr);
919 void MEDFileFieldRepresentationTree::activateTheFirst() const
921 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
922 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
923 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
925 (*it2).activateAllArrays();
930 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
932 std::size_t it0Cnt(0);
933 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
935 vtkIdType InfoOnTSId(sil->AddChild(root,edge));
936 names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
938 vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
939 std::vector<double> ts;
940 std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
941 std::size_t nbOfTS(dtits.size());
942 std::ostringstream oss3; oss3 << nbOfTS;
943 names.push_back(oss3.str());
944 for(std::size_t i=0;i<nbOfTS;i++)
946 std::ostringstream oss4; oss4 << dtits[i].first;
947 vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
948 names.push_back(oss4.str());
949 std::ostringstream oss5; oss5 << dtits[i].second;
950 vtkIdType ItId(sil->AddChild(DtId,edge));
951 names.push_back(oss5.str());
952 std::ostringstream oss6; oss6 << ts[i];
953 sil->AddChild(ItId,edge);
954 names.push_back(oss6.str());
957 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
958 std::string tsName(oss.str());
959 vtkIdType typeId0(sil->AddChild(root,edge));
960 names.push_back(tsName);
961 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
963 std::string meshName((*it1)[0].getMeshName());
964 vtkIdType typeId1(sil->AddChild(typeId0,edge));
965 names.push_back(meshName);
966 std::size_t it2Cnt(0);
967 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
969 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
970 std::string comSupStr(oss2.str());
971 vtkIdType typeId2(sil->AddChild(typeId1,edge));
972 names.push_back(comSupStr);
973 (*it2).feedSIL(_ms,meshName,sil,typeId2,edge,names);
979 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
981 int dummy0(0),dummy1(0),dummy2(0);
982 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
983 std::string ret(leaf.getMeshName());
986 for(;i<_ms->getNumberOfMeshes();i++)
988 m=_ms->getMeshAtPos(i);
989 if(m->getName()==ret)
992 if(i==_ms->getNumberOfMeshes())
993 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
994 vtkIdType typeId0(sil->AddChild(root,edge));
995 names.push_back(m->getName());
997 vtkIdType typeId1(sil->AddChild(typeId0,edge));
998 names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
999 std::vector<std::string> grps(m->getGroupsNames());
1000 for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
1002 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1003 names.push_back(*it0);
1004 std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
1005 for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
1007 sil->AddChild(typeId2,edge);
1008 names.push_back((*it1).c_str());
1012 vtkIdType typeId11(sil->AddChild(typeId0,edge));
1013 names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
1014 std::vector<std::string> fams(m->getFamiliesNames());
1015 for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
1017 sil->AddChild(typeId11,edge);
1018 int famId(m->getFamilyId((*it00).c_str()));
1019 std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
1020 names.push_back(oss.str());
1025 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
1027 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1028 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1029 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1030 if((*it2).containId(id))
1031 return (*it2).getLeafArr(id);
1032 throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
1035 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
1037 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1038 return elt.getZeName();
1041 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1043 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1044 return elt.getStatus();
1047 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1050 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1051 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1052 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1053 if((*it2).containZeName(name,ret))
1055 std::ostringstream msg; msg << "MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name \"" << name << "\" !";
1056 throw INTERP_KERNEL::Exception(msg.str().c_str());
1059 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1061 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1062 bool ret(elt.setStatus(status));//to be implemented
1066 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1069 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1070 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1071 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1072 ret=std::max(ret,(*it2).getNumberOfTS());
1079 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv)
1083 _ms=MEDFileMeshes::New(fileName);
1084 _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
1088 MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::SauvReader> sr(ParaMEDMEM::SauvReader::New(fileName));
1089 MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileData> mfd(sr->loadInMEDFileDS());
1090 _ms=mfd->getMeshes(); _ms->incrRef();
1091 int nbMeshes(_ms->getNumberOfMeshes());
1092 for(int i=0;i<nbMeshes;i++)
1094 ParaMEDMEM::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1095 ParaMEDMEM::MEDFileUMesh *tmp2(dynamic_cast<ParaMEDMEM::MEDFileUMesh *>(tmp));
1097 tmp2->forceComputationOfParts();
1099 _fields=mfd->getFields();
1100 if((ParaMEDMEM::MEDFileFields *)_fields)
1103 if(!((ParaMEDMEM::MEDFileFields *)_fields))
1105 _fields=BuildFieldFromMeshes(_ms);
1109 AppendFieldFromMeshes(_ms,_fields);
1111 _fields->removeFieldsWithoutAnyTimeStep();
1112 std::vector<std::string> meshNames(_ms->getMeshesNames());
1113 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> > fields_per_mesh(meshNames.size());
1114 for(std::size_t i=0;i<meshNames.size();i++)
1116 fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1118 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1120 for(std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1122 for(int j=0;j<(*fields)->getNumberOfFields();j++)
1124 MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1125 std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1126 allFMTSLeavesToDisplaySafe.insert(allFMTSLeavesToDisplaySafe.end(),tmp.begin(),tmp.end());
1129 std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1130 for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1132 allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1134 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1135 // memory safety part
1136 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1137 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1139 allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1140 for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1142 allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1143 allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1146 // end of memory safety part
1147 // 1st : timesteps, 2nd : meshName, 3rd : common support
1148 this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1149 for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1151 vtksys_stl::vector< vtksys_stl::string > meshNamesLoc;
1152 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1153 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1155 std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1156 vtksys_stl::vector< vtksys_stl::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1157 if(it==meshNamesLoc.end())
1159 meshNamesLoc.push_back(meshName);
1160 splitByMeshName.resize(splitByMeshName.size()+1);
1161 splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1164 splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1166 _data_structure[i].resize(meshNamesLoc.size());
1167 for(std::size_t j=0;j<splitByMeshName.size();j++)
1169 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFastCellSupportComparator> > fsp;
1170 std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1171 for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1172 sbmn[k]=splitByMeshName[j][k];
1173 //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1174 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1175 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1176 this->_data_structure[i][j].resize(commonSupSplit.size());
1177 for(std::size_t k=0;k<commonSupSplit.size();k++)
1179 commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1180 for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1182 commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1183 commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1186 for(std::size_t k=0;k<commonSupSplit.size();k++)
1187 this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1190 this->removeEmptyLeaves();
1192 this->computeFullNameInLeaves();
1195 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1197 std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1198 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1200 std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1201 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1203 std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1204 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1206 newSD1.push_back(*it2);
1208 newSD0.push_back(newSD1);
1211 newSD.push_back(newSD0);
1215 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1217 if(compInfos.size()!=1)
1219 return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1222 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1224 return _data_structure[0][0][0].getMeshName();
1227 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1230 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1231 return leaf.getTimeSteps(tk);
1234 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const
1237 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1238 meshName=leaf.getMeshName();
1239 std::vector<double> ts(leaf.getTimeSteps(tk));
1240 std::size_t zeTimeId(0);
1243 std::vector<double> ts2(ts.size());
1244 std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1245 std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1246 zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1249 if(zeTimeId==(int)ts.size())
1250 zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1251 if(zeTimeId==(int)ts.size())
1252 {//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.
1253 //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1255 double valAttachedToPos(-std::numeric_limits<double>::max());
1256 for(std::size_t i=0;i<ts.size();i++)
1260 if(ts[i]>valAttachedToPos)
1263 valAttachedToPos=ts[i];
1268 {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1269 valAttachedToPos=std::numeric_limits<double>::max();
1270 for(std::size_t i=0;i<ts.size();i++)
1272 if(ts[i]<valAttachedToPos)
1275 valAttachedToPos=ts[i];
1280 std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1281 std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1282 oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1283 std::cerr << oss.str() << std::endl;
1287 tr=new MEDStdTimeReq((int)zeTimeId);
1289 tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
1290 vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms));
1295 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1297 int nbOfActivated(0);
1298 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1299 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1300 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1301 if((*it2).isActivated())
1303 if(nbOfActivated!=1)
1305 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1306 throw INTERP_KERNEL::Exception(oss.str().c_str());
1308 int i0(0),i1(0),i2(0);
1309 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1310 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1311 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1312 if((*it2).isActivated())
1314 lev0=i0; lev1=i1; lev2=i2;
1317 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1320 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1323 os << "#############################################" << std::endl;
1324 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1327 os << "TS" << i << std::endl;
1328 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1331 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1334 os << " " << (*it2).getMeshName() << std::endl;
1335 os << " Comp" << k << std::endl;
1336 (*it2).printMySelf(os);
1340 os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1343 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms, ParaMEDMEM::MEDFileFields *ret)
1345 for(int i=0;i<ms->getNumberOfMeshes();i++)
1347 MEDFileMesh *mm(ms->getMeshAtPos(i));
1348 std::vector<int> levs(mm->getNonEmptyLevels());
1349 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileField1TS> f1tsMultiLev(ParaMEDMEM::MEDFileField1TS::New());
1350 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1353 for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1355 std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1356 for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1358 ParaMEDMEM::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1359 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1361 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1362 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1365 f->setName(mm->getName());
1366 f1tsMultiLev->setFieldNoProfileSBT(f);
1371 std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1372 if(levsExt.size()==levs.size()+1)
1374 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(1));
1375 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_NODES));
1377 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1378 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1379 arr->iota(); f->setArray(arr);
1380 f->setName(mm->getName());
1381 f1tsMultiLev->setFieldNoProfileSBT(f);
1389 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(0));
1390 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1392 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1393 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1396 f->setName(mm->getName());
1397 f1tsMultiLev->setFieldNoProfileSBT(f);
1400 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFieldMultiTS> fmtsMultiLev(ParaMEDMEM::MEDFileFieldMultiTS::New());
1401 fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1402 ret->pushField(fmtsMultiLev);
1406 ParaMEDMEM::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms)
1408 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFields> ret(ParaMEDMEM::MEDFileFields::New());
1409 AppendFieldFromMeshes(ms,ret);
1413 std::vector<std::string> MEDFileFieldRepresentationTree::SplitFieldNameIntoParts(const std::string& fullFieldName, char sep)
1415 std::vector<std::string> ret;
1417 while(pos!=std::string::npos)
1419 std::size_t curPos(fullFieldName.find_first_of(sep,pos));
1420 std::string elt(fullFieldName.substr(pos,curPos!=std::string::npos?curPos-pos:std::string::npos));
1422 pos=fullFieldName.find_first_not_of(sep,curPos);
1428 * Here the non regression tests.
1429 * const char inp0[]="";
1430 * const char exp0[]="";
1431 * const char inp1[]="field";
1432 * const char exp1[]="field";
1433 * const char inp2[]="_________";
1434 * const char exp2[]="_________";
1435 * const char inp3[]="field_p";
1436 * const char exp3[]="field_p";
1437 * const char inp4[]="field__p";
1438 * const char exp4[]="field_p";
1439 * const char inp5[]="field_p__";
1440 * const char exp5[]="field_p";
1441 * const char inp6[]="field_p_";
1442 * const char exp6[]="field_p";
1443 * const char inp7[]="field_____EDFGEG//sdkjf_____PP_______________";
1444 * const char exp7[]="field_EDFGEG//sdkjf_PP";
1445 * const char inp8[]="field_____EDFGEG//sdkjf_____PP";
1446 * const char exp8[]="field_EDFGEG//sdkjf_PP";
1447 * const char inp9[]="_field_____EDFGEG//sdkjf_____PP_______________";
1448 * const char exp9[]="field_EDFGEG//sdkjf_PP";
1449 * const char inp10[]="___field_____EDFGEG//sdkjf_____PP_______________";
1450 * const char exp10[]="field_EDFGEG//sdkjf_PP";
1452 std::string MEDFileFieldRepresentationTree::PostProcessFieldName(const std::string& fullFieldName)
1454 static const char SEP('_');
1455 std::vector<std::string> v(SplitFieldNameIntoParts(fullFieldName,SEP));
1457 return fullFieldName;//should never happen
1461 return fullFieldName;
1465 std::string ret(v[0]);
1466 for(std::size_t i=1;i<v.size();i++)
1471 { ret+=SEP; ret+=v[i]; }
1477 return fullFieldName;
1483 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1487 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1492 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1494 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1496 throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1502 * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1503 * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1504 * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1506 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1508 std::size_t sz(ts.size());
1509 bool isInHumanRange(true);
1511 for(std::size_t i=0;i<sz;i++)
1514 if(ts[i]<=-1e299 || ts[i]>=1e299)
1515 isInHumanRange=false;
1518 return processedUsingPairOfIds(tsPairs);
1520 return processedUsingPairOfIds(tsPairs);
1521 _postprocessed_time=ts;
1522 return getPostProcessedTime();
1527 * idem than 0, except that ts is preaccumulated before invoking policy 0.
1529 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1531 std::size_t sz(ts.size());
1532 std::vector<double> ts2(sz);
1534 for(std::size_t i=0;i<sz;i++)
1539 return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1542 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1544 std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1548 void TimeKeeper::printSelf(std::ostream& oss) const
1550 std::size_t sz(_activated_ts.size());
1551 for(std::size_t i=0;i<sz;i++)
1553 oss << "(" << i << "," << _activated_ts[i].first << "), ";
1557 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1559 std::size_t sz(_activated_ts.size());
1560 std::vector<bool> ret(sz);
1561 for(std::size_t i=0;i<sz;i++)
1563 ret[i]=_activated_ts[i].first;
1568 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1570 std::size_t sz(tsPairs.size());
1571 std::set<int> s0,s1;
1572 for(std::size_t i=0;i<sz;i++)
1573 { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1576 _postprocessed_time.resize(sz);
1577 for(std::size_t i=0;i<sz;i++)
1578 _postprocessed_time[i]=(double)tsPairs[i].first;
1579 return getPostProcessedTime();
1583 _postprocessed_time.resize(sz);
1584 for(std::size_t i=0;i<sz;i++)
1585 _postprocessed_time[i]=(double)tsPairs[i].second;
1586 return getPostProcessedTime();
1588 //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1589 _postprocessed_time.resize(sz);
1590 for(std::size_t i=0;i<sz;i++)
1591 _postprocessed_time[i]=(double)i;
1592 return getPostProcessedTime();
1595 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1597 _activated_ts.resize(maxNumberOfTS);
1598 for(int i=0;i<maxNumberOfTS;i++)
1600 std::ostringstream oss; oss << "000" << i;
1601 _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());