1 // Copyright (C) 2010-2015 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 #ifdef MEDREADER_USE_MPI
32 #include "ParaMEDFileMesh.hxx"
35 #include "vtkXMLUnstructuredGridWriter.h"//
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"
53 #include "vtksys/stl/string"
54 #include "vtksys/ios/fstream"
55 #include "vtksys/stl/algorithm"
56 #include "vtkMutableDirectedGraph.h"
58 using namespace ParaMEDMEM;
60 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
62 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
64 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
66 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
68 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
70 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
72 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
74 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
76 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
78 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
80 vtkIdTypeArray *ELGACmp::findOrCreate(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, bool& isNew) const
82 vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
91 return createNew(globs,locsReallyUsed,vtkd,ds);
95 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd) const
97 std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
98 if(it==_loc_names.end())
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++)
105 key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
107 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
111 vtkIdTypeArray *ELGACmp::createNew(const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds) const
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;
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++)
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++)
142 const double *pt0(calculator.getFunctionValues(i));
143 if(ct!=INTERP_KERNEL::NORM_HEXA27)
144 std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
147 for(int j=0;j<27;j++)
148 shape[nbPtsPerCell*i+j]=pt0[MEDMeshMultiLev::HEXA27_PERM_ARRAY[j]];
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]));
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));
160 vtkIdType ncell(ds->GetNumberOfCells());
161 int *pt(new int[ncell]),offset(0);
162 for(vtkIdType cellId=0;cellId<ncell;cellId++)
164 vtkCell *cell(ds->GetCell(cellId));
165 int delta(m[cell->GetCellType()]);
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);
180 _defs.push_back(defs);
183 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
185 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
186 ds->GetCellData()->AddArray(*it);
191 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
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();
200 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
204 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>& arr):ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
206 std::vector< std::vector<ParaMEDMEM::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
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();
216 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
218 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileAnyTypeFieldMultiTS>::operator=(other);
221 _ze_name=other._ze_name;
222 _ze_full_name.clear();
226 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
231 int MEDFileFieldRepresentationLeavesArrays::getId() const
236 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
238 return _ze_full_name;
241 const char *MEDFileFieldRepresentationLeavesArrays::getZeNameC() const
243 return _ze_full_name.c_str();
246 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
248 vtkIdType refId(sil->AddChild(root,edge));
249 names.push_back(_ze_name);
251 if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
253 sil->AddChild(refId,edge);
254 names.push_back(std::string());
258 void MEDFileFieldRepresentationLeavesArrays::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
260 std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
261 _ze_full_name=oss3.str();
264 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
269 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
271 bool ret(_activated!=status);
276 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshStruct *mst, vtkDataSet *ds) const
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());
282 for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
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);
290 crudeArr=f1tsPtrDbl->getUndergroundDataArray();
292 crudeArr=f1tsPtrInt->getUndergroundDataArray();
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));
300 std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
302 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
303 vtkFieldData *att(0);
308 att=ds->GetCellData();
313 att=ds->GetPointData();
318 att=ds->GetFieldData();
323 att=ds->GetFieldData();
327 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
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)
338 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vd->accessToMemArray().setSpecificDeallocator(0);
342 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
344 std::string name(tr->buildName(f1ts->getName()));
345 vtkd->SetName(name.c_str());
348 if(discs[0]==ON_GAUSS_PT)
351 _elga_cmp.findOrCreate(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp);
353 if(discs[0]==ON_GAUSS_NE)
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++)
362 vtkCell *cell(ds->GetCell(cellId));
363 int delta(cell->GetNumberOfPoints());
364 cellTypes.insert(cell->GetCellType());
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);
376 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
377 for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
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)
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]);
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++)
396 const double *pt0(calculator.getFunctionValues(i));
397 std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
399 def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
401 key->Set(elno->GetInformation(),def,*it);
402 key->Set(vtkd->GetInformation(),def,*it);
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)
418 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vi->accessToMemArray().setSpecificDeallocator(0);
422 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
424 std::string name(tr->buildName(f1ts->getName()));
425 vtkd->SetName(name.c_str());
430 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
434 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
436 _elga_cmp.appendELGAIfAny(ds);
441 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
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)
448 for(std::size_t i=0;i<arr.size();i++)
449 _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
452 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
455 _cached_ds->Delete();
458 bool MEDFileFieldRepresentationLeaves::empty() const
460 const MEDFileFastCellSupportComparator *fcscp(_fsp);
461 return fcscp==0 || _arrays.empty();
464 void MEDFileFieldRepresentationLeaves::setId(int& id) const
466 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
470 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
472 return _arrays[0]->getMeshName();
475 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
477 return (int)_arrays.size();
480 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
482 return _arrays[0]->getNumberOfTS();
485 void MEDFileFieldRepresentationLeaves::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
487 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
488 (*it).computeFullNameInLeaves(tsName,meshName,comSupStr);
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.
494 void MEDFileFieldRepresentationLeaves::feedSIL(const ParaMEDMEM::MEDFileMeshes *ms, const std::string& meshName, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
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);
501 vtkIdType root3(sil->AddChild(root,edge));
502 names.push_back(std::string("InfoOnGeoType"));
503 const ParaMEDMEM::MEDFileMesh *m(0);
505 m=ms->getMeshWithName(meshName);
506 const ParaMEDMEM::MEDFileFastCellSupportComparator *fsp(_fsp);
507 if(!fsp || fsp->getNumberOfTS()==0)
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++)
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);
519 bool MEDFileFieldRepresentationLeaves::containId(int id) const
521 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
522 if((*it).getId()==id)
527 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
529 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
530 if((*it).getZeName()==name)
538 void MEDFileFieldRepresentationLeaves::dumpState(std::map<std::string,bool>& status) const
540 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
541 status[(*it).getZeName()]=(*it).getStatus();
544 bool MEDFileFieldRepresentationLeaves::isActivated() const
546 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
547 if((*it).getStatus())
552 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
554 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
556 os << " - " << (*it0).getZeName() << " (";
557 if((*it0).getStatus())
561 os << ")" << std::endl;
565 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
567 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
568 (*it).setStatus(true);
571 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
573 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
574 if((*it).getId()==id)
576 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
579 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
582 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
583 std::vector<double> ret;
584 std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
585 return tk.getTimeStepsRegardingPolicy(dtits,ret);
588 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
591 return _arrays[0]->getTimeSteps(ts);
595 return std::vector< std::pair<int,int> >();
599 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
601 std::ostringstream oss;
602 oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
603 std::vector<double> ret1;
604 std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
605 std::size_t sz(ret1.size());
606 for(std::size_t i=0;i<sz;i++)
608 oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
611 std::string tmp(oss.str());
612 if(tmp.size()>200 && i!=sz-1)
622 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const ParaMEDMEM::MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDMeshMultiLev *mml, const ParaMEDMEM::MEDFileMeshes *meshes, vtkDataSet *ds) const
625 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
626 MEDCouplingAutoRefCountObjectPtr<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
627 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
628 if((*it).getStatus())
630 (*it).appendFields(tr,globs,mml,mst,ds);
631 (*it).appendELGAIfAny(ds);
635 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
637 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
638 DataArrayDouble *coordsMC(0);
639 DataArrayByte *typesMC(0);
640 DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
641 bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
642 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coordsSafe(coordsMC);
643 MEDCouplingAutoRefCountObjectPtr<DataArrayByte> typesSafe(typesMC);
644 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
646 int nbOfCells(typesSafe->getNbOfElems());
647 vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
648 vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
649 cellTypes->SetArray(reinterpret_cast<unsigned char *>(typesSafe->getPointer()),nbOfCells,0,VTK_DATA_ARRAY_FREE); typesSafe->accessToMemArray().setSpecificDeallocator(0);
650 vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
651 cellLocations->SetArray(cellLocationsSafe->getPointer(),nbOfCells,0,VTK_DATA_ARRAY_FREE); cellLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
652 vtkCellArray *cells(vtkCellArray::New());
653 vtkIdTypeArray *cells2(vtkIdTypeArray::New());
654 cells2->SetArray(cellsSafe->getPointer(),cellsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); cellsSafe->accessToMemArray().setSpecificDeallocator(0);
655 cells->SetCells(nbOfCells,cells2);
657 if(faceLocationsMC!=0 && facesMC!=0)
659 vtkIdTypeArray *faces(vtkIdTypeArray::New());
660 faces->SetArray(facesSafe->getPointer(),facesSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); facesSafe->accessToMemArray().setSpecificDeallocator(0);
661 vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
662 faceLocations->SetArray(faceLocationsSafe->getPointer(),faceLocationsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); faceLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
663 ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
664 faceLocations->Delete();
668 ret->SetCells(cellTypes,cellLocations,cells);
670 cellLocations->Delete();
672 vtkPoints *pts(vtkPoints::New());
673 vtkDoubleArray *pts2(vtkDoubleArray::New());
674 pts2->SetNumberOfComponents(3);
677 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);
678 coordsSafe->accessToMemArray().setSpecificDeallocator(0);
681 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
690 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(ParaMEDMEM::MEDCMeshMultiLev *mm) const
692 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
694 std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
695 vtkDoubleArray *vtkTmp(0);
696 vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
697 std::size_t dim(arrs.size());
699 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
700 int sizePerAxe[3]={1,1,1};
701 sizePerAxe[0]=arrs[0]->getNbOfElems();
703 sizePerAxe[1]=arrs[1]->getNbOfElems();
705 sizePerAxe[2]=arrs[2]->getNbOfElems();
706 ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
707 vtkTmp=vtkDoubleArray::New();
708 vtkTmp->SetNumberOfComponents(1);
710 vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
712 { vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[0]->accessToMemArray().setSpecificDeallocator(0); }
713 ret->SetXCoordinates(vtkTmp);
718 vtkTmp=vtkDoubleArray::New();
719 vtkTmp->SetNumberOfComponents(1);
721 vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
723 { vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[1]->accessToMemArray().setSpecificDeallocator(0); }
724 ret->SetYCoordinates(vtkTmp);
730 vtkTmp=vtkDoubleArray::New();
731 vtkTmp->SetNumberOfComponents(1);
733 vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
735 { vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[2]->accessToMemArray().setSpecificDeallocator(0); }
736 ret->SetZCoordinates(vtkTmp);
743 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(ParaMEDMEM::MEDCurveLinearMeshMultiLev *mm) const
745 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
746 int meshStr[3]={1,1,1};
747 DataArrayDouble *coords(0);
748 std::vector<int> nodeStrct;
750 mm->buildVTUArrays(coords,nodeStrct,isInternal);
751 std::size_t dim(nodeStrct.size());
753 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
754 meshStr[0]=nodeStrct[0];
756 meshStr[1]=nodeStrct[1];
758 meshStr[2]=nodeStrct[2];
759 vtkStructuredGrid *ret(vtkStructuredGrid::New());
760 ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
761 vtkDoubleArray *da(vtkDoubleArray::New());
762 da->SetNumberOfComponents(3);
763 if(coords->getNumberOfComponents()==3)
766 da->SetArray(coords->getPointer(),coords->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);//VTK has not the ownership of double * because MEDLoader main struct has it !
768 { da->SetArray(coords->getPointer(),coords->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); coords->accessToMemArray().setSpecificDeallocator(0); }
772 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
773 da->SetArray(coords2->getPointer(),coords2->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);//let VTK deal with double *
774 coords2->accessToMemArray().setSpecificDeallocator(0);
777 vtkPoints *points=vtkPoints::New();
778 ret->SetPoints(points);
785 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const ParaMEDMEM::MEDFileMeshes *meshes) const
787 static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
789 //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
790 MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change
791 MEDCouplingAutoRefCountObjectPtr<MEDMeshMultiLev> mml2(mml->prepare());
792 MEDMeshMultiLev *ptMML2(mml2);
795 MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
796 MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
797 MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
801 ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
805 ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
809 ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
812 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
813 _cached_ds=ret->NewInstance();
814 _cached_ds->ShallowCopy(ret);
818 ret=_cached_ds->NewInstance();
819 ret->ShallowCopy(_cached_ds);
822 appendFields(tr,globs,mml,meshes,ret);
823 // The arrays links to mesh
824 DataArrayInt *famCells(0),*numCells(0);
825 bool noCpyFamCells(false),noCpyNumCells(false);
826 ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
829 vtkIntArray *vtkTab(vtkIntArray::New());
830 vtkTab->SetNumberOfComponents(1);
831 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
833 vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
835 { vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famCells->accessToMemArray().setSpecificDeallocator(0); }
836 ret->GetCellData()->AddArray(vtkTab);
840 ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
843 vtkIntArray *vtkTab(vtkIntArray::New());
844 vtkTab->SetNumberOfComponents(1);
845 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
847 vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
849 { vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numCells->accessToMemArray().setSpecificDeallocator(0); }
850 ret->GetCellData()->AddArray(vtkTab);
854 // The arrays links to mesh
855 DataArrayInt *famNodes(0),*numNodes(0);
856 bool noCpyFamNodes(false),noCpyNumNodes(false);
857 ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
860 vtkIntArray *vtkTab(vtkIntArray::New());
861 vtkTab->SetNumberOfComponents(1);
862 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
864 vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
866 { vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famNodes->accessToMemArray().setSpecificDeallocator(0); }
867 ret->GetPointData()->AddArray(vtkTab);
871 ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
874 vtkIntArray *vtkTab(vtkIntArray::New());
875 vtkTab->SetNumberOfComponents(1);
876 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
878 vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
880 { vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numNodes->accessToMemArray().setSpecificDeallocator(0); }
881 ret->GetPointData()->AddArray(vtkTab);
888 //////////////////////
890 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
894 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
897 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
898 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
899 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
900 ret+=(*it2).getNumberOfArrays();
904 void MEDFileFieldRepresentationTree::assignIds() const
907 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
908 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
909 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
913 void MEDFileFieldRepresentationTree::computeFullNameInLeaves() const
915 std::size_t it0Cnt(0);
916 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
918 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
919 std::string tsName(oss.str());
920 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
922 std::string meshName((*it1)[0].getMeshName());
923 std::size_t it2Cnt(0);
924 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
926 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
927 std::string comSupStr(oss2.str());
928 (*it2).computeFullNameInLeaves(tsName,meshName,comSupStr);
934 void MEDFileFieldRepresentationTree::activateTheFirst() const
936 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
937 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
938 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
940 (*it2).activateAllArrays();
945 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
947 std::size_t it0Cnt(0);
948 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
950 vtkIdType InfoOnTSId(sil->AddChild(root,edge));
951 names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
953 vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
954 std::vector<double> ts;
955 std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
956 std::size_t nbOfTS(dtits.size());
957 std::ostringstream oss3; oss3 << nbOfTS;
958 names.push_back(oss3.str());
959 for(std::size_t i=0;i<nbOfTS;i++)
961 std::ostringstream oss4; oss4 << dtits[i].first;
962 vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
963 names.push_back(oss4.str());
964 std::ostringstream oss5; oss5 << dtits[i].second;
965 vtkIdType ItId(sil->AddChild(DtId,edge));
966 names.push_back(oss5.str());
967 std::ostringstream oss6; oss6 << ts[i];
968 sil->AddChild(ItId,edge);
969 names.push_back(oss6.str());
972 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
973 std::string tsName(oss.str());
974 vtkIdType typeId0(sil->AddChild(root,edge));
975 names.push_back(tsName);
976 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
978 std::string meshName((*it1)[0].getMeshName());
979 vtkIdType typeId1(sil->AddChild(typeId0,edge));
980 names.push_back(meshName);
981 std::size_t it2Cnt(0);
982 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
984 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
985 std::string comSupStr(oss2.str());
986 vtkIdType typeId2(sil->AddChild(typeId1,edge));
987 names.push_back(comSupStr);
988 (*it2).feedSIL(_ms,meshName,sil,typeId2,edge,names);
994 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
996 int dummy0(0),dummy1(0),dummy2(0);
997 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
998 std::string ret(leaf.getMeshName());
1001 for(;i<_ms->getNumberOfMeshes();i++)
1003 m=_ms->getMeshAtPos(i);
1004 if(m->getName()==ret)
1007 if(i==_ms->getNumberOfMeshes())
1008 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
1009 vtkIdType typeId0(sil->AddChild(root,edge));
1010 names.push_back(m->getName());
1012 vtkIdType typeId1(sil->AddChild(typeId0,edge));
1013 names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
1014 std::vector<std::string> grps(m->getGroupsNames());
1015 for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
1017 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1018 names.push_back(*it0);
1019 std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
1020 for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
1022 sil->AddChild(typeId2,edge);
1023 names.push_back((*it1).c_str());
1027 vtkIdType typeId11(sil->AddChild(typeId0,edge));
1028 names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
1029 std::vector<std::string> fams(m->getFamiliesNames());
1030 for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
1032 sil->AddChild(typeId11,edge);
1033 int famId(m->getFamilyId((*it00).c_str()));
1034 std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
1035 names.push_back(oss.str());
1040 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
1042 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1043 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1044 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1045 if((*it2).containId(id))
1046 return (*it2).getLeafArr(id);
1047 throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
1050 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
1052 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1053 return elt.getZeName();
1056 const char *MEDFileFieldRepresentationTree::getNameOfC(int id) const
1058 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1059 return elt.getZeNameC();
1062 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1064 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1065 return elt.getStatus();
1068 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1071 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1072 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1073 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1074 if((*it2).containZeName(name,ret))
1076 std::ostringstream msg; msg << "MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name \"" << name << "\" !";
1077 throw INTERP_KERNEL::Exception(msg.str().c_str());
1080 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1082 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1083 bool ret(elt.setStatus(status));//to be implemented
1087 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1090 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1091 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1092 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1093 ret=std::max(ret,(*it2).getNumberOfTS());
1100 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts)
1104 if((iPart==-1 && nbOfParts==-1) || (iPart==0 && nbOfParts==1))
1106 _ms=MEDFileMeshes::New(fileName);
1107 _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
1111 #ifdef MEDREADER_USE_MPI
1112 _ms=ParaMEDFileMeshes::New(iPart,nbOfParts,fileName);
1113 int nbMeshes(_ms->getNumberOfMeshes());
1114 for(int i=0;i<nbMeshes;i++)
1116 ParaMEDMEM::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1117 ParaMEDMEM::MEDFileUMesh *tmp2(dynamic_cast<ParaMEDMEM::MEDFileUMesh *>(tmp));
1119 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp3(tmp2->zipCoords());
1121 _fields=MEDFileFields::LoadPartOf(fileName,false,_ms);//false is important to not read the values
1123 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::loadMainStructureOfFile : request for iPart/nbOfParts=" << iPart << "/" << nbOfParts << " whereas Plugin not compiled with MPI !";
1124 throw INTERP_KERNEL::Exception(oss.str().c_str());
1130 MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::SauvReader> sr(ParaMEDMEM::SauvReader::New(fileName));
1131 MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileData> mfd(sr->loadInMEDFileDS());
1132 _ms=mfd->getMeshes(); _ms->incrRef();
1133 int nbMeshes(_ms->getNumberOfMeshes());
1134 for(int i=0;i<nbMeshes;i++)
1136 ParaMEDMEM::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1137 ParaMEDMEM::MEDFileUMesh *tmp2(dynamic_cast<ParaMEDMEM::MEDFileUMesh *>(tmp));
1139 tmp2->forceComputationOfParts();
1141 _fields=mfd->getFields();
1142 if((ParaMEDMEM::MEDFileFields *)_fields)
1145 if(!((ParaMEDMEM::MEDFileFields *)_fields))
1147 _fields=BuildFieldFromMeshes(_ms);
1151 AppendFieldFromMeshes(_ms,_fields);
1153 _fields->removeFieldsWithoutAnyTimeStep();
1154 std::vector<std::string> meshNames(_ms->getMeshesNames());
1155 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> > fields_per_mesh(meshNames.size());
1156 for(std::size_t i=0;i<meshNames.size();i++)
1158 fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1160 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1162 for(std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1164 for(int j=0;j<(*fields)->getNumberOfFields();j++)
1166 MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1167 std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1169 for(std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > >::const_iterator it=tmp.begin();it!=tmp.end();it++)
1171 if(!(*it)->presenceOfMultiDiscPerGeoType())
1172 allFMTSLeavesToDisplaySafe.push_back(*it);
1174 {// The case of some parts of field have more than one discretization per geo type.
1175 std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > > subTmp((*it)->splitMultiDiscrPerGeoTypes());
1176 std::size_t it0Cnt(0);
1177 for(std::vector< MEDCouplingAutoRefCountObjectPtr< MEDFileAnyTypeFieldMultiTS > >::iterator it0=subTmp.begin();it0!=subTmp.end();it0++,it0Cnt++)//not const because setName
1179 std::ostringstream oss; oss << (*it0)->getName() << "_" << std::setfill('M') << std::setw(3) << it0Cnt;
1180 (*it0)->setName(oss.str());
1181 allFMTSLeavesToDisplaySafe.push_back(*it0);
1188 std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1189 for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1191 allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1193 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1194 // memory safety part
1195 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1196 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1198 allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1199 for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1201 allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1202 allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1205 // end of memory safety part
1206 // 1st : timesteps, 2nd : meshName, 3rd : common support
1207 this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1208 for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1210 vtksys_stl::vector< vtksys_stl::string > meshNamesLoc;
1211 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1212 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1214 std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1215 vtksys_stl::vector< vtksys_stl::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1216 if(it==meshNamesLoc.end())
1218 meshNamesLoc.push_back(meshName);
1219 splitByMeshName.resize(splitByMeshName.size()+1);
1220 splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1223 splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1225 _data_structure[i].resize(meshNamesLoc.size());
1226 for(std::size_t j=0;j<splitByMeshName.size();j++)
1228 std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileFastCellSupportComparator> > fsp;
1229 std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1230 for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1231 sbmn[k]=splitByMeshName[j][k];
1232 //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1233 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1234 std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1235 this->_data_structure[i][j].resize(commonSupSplit.size());
1236 for(std::size_t k=0;k<commonSupSplit.size();k++)
1238 commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1239 for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1241 commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1242 commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1245 for(std::size_t k=0;k<commonSupSplit.size();k++)
1246 this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1249 this->removeEmptyLeaves();
1251 this->computeFullNameInLeaves();
1254 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1256 std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1257 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1259 std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1260 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1262 std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1263 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1265 newSD1.push_back(*it2);
1267 newSD0.push_back(newSD1);
1270 newSD.push_back(newSD0);
1274 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1276 if(compInfos.size()!=1)
1278 return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1281 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1283 return _data_structure[0][0][0].getMeshName();
1286 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1289 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1290 return leaf.getTimeSteps(tk);
1293 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const
1296 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1297 meshName=leaf.getMeshName();
1298 std::vector<double> ts(leaf.getTimeSteps(tk));
1299 std::size_t zeTimeId(0);
1302 std::vector<double> ts2(ts.size());
1303 std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1304 std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1305 zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1308 if(zeTimeId==(int)ts.size())
1309 zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1310 if(zeTimeId==(int)ts.size())
1311 {//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.
1312 //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1314 double valAttachedToPos(-std::numeric_limits<double>::max());
1315 for(std::size_t i=0;i<ts.size();i++)
1319 if(ts[i]>valAttachedToPos)
1322 valAttachedToPos=ts[i];
1327 {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1328 valAttachedToPos=std::numeric_limits<double>::max();
1329 for(std::size_t i=0;i<ts.size();i++)
1331 if(ts[i]<valAttachedToPos)
1334 valAttachedToPos=ts[i];
1339 std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1340 std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1341 oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1342 std::cerr << oss.str() << std::endl;
1346 tr=new MEDStdTimeReq((int)zeTimeId);
1348 tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
1349 vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms));
1354 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1356 int nbOfActivated(0);
1357 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1358 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1359 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1360 if((*it2).isActivated())
1362 if(nbOfActivated!=1)
1364 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1365 throw INTERP_KERNEL::Exception(oss.str().c_str());
1367 int i0(0),i1(0),i2(0);
1368 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1369 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1370 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1371 if((*it2).isActivated())
1373 lev0=i0; lev1=i1; lev2=i2;
1376 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1379 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1382 os << "#############################################" << std::endl;
1383 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1386 os << "TS" << i << std::endl;
1387 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1390 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1393 os << " " << (*it2).getMeshName() << std::endl;
1394 os << " Comp" << k << std::endl;
1395 (*it2).printMySelf(os);
1399 os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1402 std::map<std::string,bool> MEDFileFieldRepresentationTree::dumpState() const
1404 std::map<std::string,bool> ret;
1405 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1406 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1407 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1408 (*it2).dumpState(ret);
1412 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms, ParaMEDMEM::MEDFileFields *ret)
1415 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::AppendFieldFromMeshes : internal error ! NULL ret !");
1416 for(int i=0;i<ms->getNumberOfMeshes();i++)
1418 MEDFileMesh *mm(ms->getMeshAtPos(i));
1419 std::vector<int> levs(mm->getNonEmptyLevels());
1420 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileField1TS> f1tsMultiLev(ParaMEDMEM::MEDFileField1TS::New());
1421 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1424 for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1426 std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1427 for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1429 ParaMEDMEM::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1430 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1432 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1433 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1436 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1437 f1tsMultiLev->setFieldNoProfileSBT(f);
1442 std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1443 if(levsExt.size()==levs.size()+1)
1445 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(1));
1446 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_NODES));
1448 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1449 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1450 arr->iota(); f->setArray(arr);
1451 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1452 f1tsMultiLev->setFieldNoProfileSBT(f);
1460 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingMesh> m(mm->getGenMeshAtLevel(0));
1461 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDCouplingFieldDouble> f(ParaMEDMEM::MEDCouplingFieldDouble::New(ParaMEDMEM::ON_CELLS));
1463 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::DataArrayDouble> arr(ParaMEDMEM::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1464 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1467 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1468 f1tsMultiLev->setFieldNoProfileSBT(f);
1471 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFieldMultiTS> fmtsMultiLev(ParaMEDMEM::MEDFileFieldMultiTS::New());
1472 fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1473 ret->pushField(fmtsMultiLev);
1477 std::string MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh(const std::string& meshName, const ParaMEDMEM::MEDFileFields *ret)
1479 static const char KEY_STR_TO_AVOID_COLLIDE[]="MESH@";
1481 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh : internal error ! NULL ret !");
1482 std::vector<std::string> fieldNamesAlreadyExisting(ret->getFieldsNames());
1483 if(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),meshName)==fieldNamesAlreadyExisting.end())
1485 std::string tmpName(KEY_STR_TO_AVOID_COLLIDE); tmpName+=meshName;
1486 while(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),tmpName)!=fieldNamesAlreadyExisting.end())
1487 tmpName=std::string(KEY_STR_TO_AVOID_COLLIDE)+tmpName;
1491 ParaMEDMEM::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const ParaMEDMEM::MEDFileMeshes *ms)
1493 ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr<ParaMEDMEM::MEDFileFields> ret(ParaMEDMEM::MEDFileFields::New());
1494 AppendFieldFromMeshes(ms,ret);
1498 std::vector<std::string> MEDFileFieldRepresentationTree::SplitFieldNameIntoParts(const std::string& fullFieldName, char sep)
1500 std::vector<std::string> ret;
1502 while(pos!=std::string::npos)
1504 std::size_t curPos(fullFieldName.find_first_of(sep,pos));
1505 std::string elt(fullFieldName.substr(pos,curPos!=std::string::npos?curPos-pos:std::string::npos));
1507 pos=fullFieldName.find_first_not_of(sep,curPos);
1513 * Here the non regression tests.
1514 * const char inp0[]="";
1515 * const char exp0[]="";
1516 * const char inp1[]="field";
1517 * const char exp1[]="field";
1518 * const char inp2[]="_________";
1519 * const char exp2[]="_________";
1520 * const char inp3[]="field_p";
1521 * const char exp3[]="field_p";
1522 * const char inp4[]="field__p";
1523 * const char exp4[]="field_p";
1524 * const char inp5[]="field_p__";
1525 * const char exp5[]="field_p";
1526 * const char inp6[]="field_p_";
1527 * const char exp6[]="field_p";
1528 * const char inp7[]="field_____EDFGEG//sdkjf_____PP_______________";
1529 * const char exp7[]="field_EDFGEG//sdkjf_PP";
1530 * const char inp8[]="field_____EDFGEG//sdkjf_____PP";
1531 * const char exp8[]="field_EDFGEG//sdkjf_PP";
1532 * const char inp9[]="_field_____EDFGEG//sdkjf_____PP_______________";
1533 * const char exp9[]="field_EDFGEG//sdkjf_PP";
1534 * const char inp10[]="___field_____EDFGEG//sdkjf_____PP_______________";
1535 * const char exp10[]="field_EDFGEG//sdkjf_PP";
1537 std::string MEDFileFieldRepresentationTree::PostProcessFieldName(const std::string& fullFieldName)
1539 static const char SEP('_');
1540 std::vector<std::string> v(SplitFieldNameIntoParts(fullFieldName,SEP));
1542 return fullFieldName;//should never happen
1546 return fullFieldName;
1550 std::string ret(v[0]);
1551 for(std::size_t i=1;i<v.size();i++)
1556 { ret+=SEP; ret+=v[i]; }
1562 return fullFieldName;
1568 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1572 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1577 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1579 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1581 throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1587 * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1588 * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1589 * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1591 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1593 std::size_t sz(ts.size());
1594 bool isInHumanRange(true);
1596 for(std::size_t i=0;i<sz;i++)
1599 if(ts[i]<=-1e299 || ts[i]>=1e299)
1600 isInHumanRange=false;
1603 return processedUsingPairOfIds(tsPairs);
1605 return processedUsingPairOfIds(tsPairs);
1606 _postprocessed_time=ts;
1607 return getPostProcessedTime();
1612 * idem than 0, except that ts is preaccumulated before invoking policy 0.
1614 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1616 std::size_t sz(ts.size());
1617 std::vector<double> ts2(sz);
1619 for(std::size_t i=0;i<sz;i++)
1624 return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1627 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1629 std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1633 void TimeKeeper::printSelf(std::ostream& oss) const
1635 std::size_t sz(_activated_ts.size());
1636 for(std::size_t i=0;i<sz;i++)
1638 oss << "(" << i << "," << _activated_ts[i].first << "), ";
1642 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1644 std::size_t sz(_activated_ts.size());
1645 std::vector<bool> ret(sz);
1646 for(std::size_t i=0;i<sz;i++)
1648 ret[i]=_activated_ts[i].first;
1653 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1655 std::size_t sz(tsPairs.size());
1656 std::set<int> s0,s1;
1657 for(std::size_t i=0;i<sz;i++)
1658 { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1661 _postprocessed_time.resize(sz);
1662 for(std::size_t i=0;i<sz;i++)
1663 _postprocessed_time[i]=(double)tsPairs[i].first;
1664 return getPostProcessedTime();
1668 _postprocessed_time.resize(sz);
1669 for(std::size_t i=0;i<sz;i++)
1670 _postprocessed_time[i]=(double)tsPairs[i].second;
1671 return getPostProcessedTime();
1673 //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1674 _postprocessed_time.resize(sz);
1675 for(std::size_t i=0;i<sz;i++)
1676 _postprocessed_time[i]=(double)i;
1677 return getPostProcessedTime();
1680 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1682 _activated_ts.resize(maxNumberOfTS);
1683 for(int i=0;i<maxNumberOfTS;i++)
1685 std::ostringstream oss; oss << "000" << i;
1686 _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());