1 // Copyright (C) 2010-2016 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 "MEDFileBlowStrEltUp.hxx"
29 #include "MEDFileData.hxx"
30 #include "SauvReader.hxx"
32 #ifdef MEDREADER_USE_MPI
33 #include "ParaMEDFileMesh.hxx"
36 #include "vtkXMLUnstructuredGridWriter.h"//
38 #include "vtkUnstructuredGrid.h"
39 #include "vtkRectilinearGrid.h"
40 #include "vtkStructuredGrid.h"
41 #include "vtkUnsignedCharArray.h"
42 #include "vtkQuadratureSchemeDefinition.h"
43 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
44 #include "vtkInformationIntegerKey.h"
45 #include "vtkInformation.h"
46 #include "vtkDataArrayTemplate.h"
47 #include "vtkIdTypeArray.h"
48 #include "vtkDoubleArray.h"
49 #include "vtkIntArray.h"
50 #include "vtkCellArray.h"
51 #include "vtkPointData.h"
52 #include "vtkFieldData.h"
53 #include "vtkCellData.h"
55 #include "vtkMutableDirectedGraph.h"
57 using namespace MEDCoupling;
59 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
61 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
63 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
65 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
67 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
69 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
71 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
73 const char MEDFileFieldRepresentationLeavesArrays::GLOBAL_NODE_ID_NAME[]="GlobalNodeIds";// WARNING DO NOT CHANGE IT BEFORE HAVING CHECKED IN PV SOURCES !
75 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
77 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
79 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
81 vtkIdTypeArray *ELGACmp::findOrCreate(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, bool& isNew, ExportedTinyInfo *internalInfo) const
83 vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
92 return createNew(globs,locsReallyUsed,vtkd,ds,internalInfo);
96 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd) const
98 std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
99 if(it==_loc_names.end())
101 std::size_t pos(std::distance(_loc_names.begin(),it));
102 vtkIdTypeArray *ret(_elgas[pos]);
103 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
104 for(std::vector<std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it=_defs[pos].begin();it!=_defs[pos].end();it++)
106 key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
108 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
112 vtkIdTypeArray *ELGACmp::createNew(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
114 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
115 std::vector< std::vector<std::string> > locNames(_loc_names);
116 std::vector<vtkIdTypeArray *> elgas(_elgas);
117 std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > defs;
119 std::vector< std::vector<std::string> >::const_iterator it(std::find(locNames.begin(),locNames.end(),locsReallyUsed));
120 if(it!=locNames.end())
121 throw INTERP_KERNEL::Exception("ELGACmp::createNew : Method is expected to be called after isExisting call ! Entry already exists !");
122 locNames.push_back(locsReallyUsed);
123 vtkIdTypeArray *elga(vtkIdTypeArray::New());
124 elga->SetNumberOfComponents(1);
125 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
126 std::map<unsigned char,int> m;
127 for(std::vector<std::string>::const_iterator it=locsReallyUsed.begin();it!=locsReallyUsed.end();it++)
129 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
130 const MEDFileFieldLoc& loc(globs->getLocalization((*it).c_str()));
131 INTERP_KERNEL::NormalizedCellType ct(loc.getGeoType());
132 unsigned char vtkType(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE[ct]);
133 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
134 int nbGaussPt(loc.getNbOfGaussPtPerCell()),nbPtsPerCell((int)cm.getNumberOfNodes()),dimLoc(loc.getDimension());
135 // 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.
136 std::vector<double> gsCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getGaussCoords()));
137 std::vector<double> refCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getRefCoords()));
139 internalInfo->pushGaussAdditionnalInfo(vtkType,dimLoc,refCoods2,gsCoods2);
140 double *shape(new double[nbPtsPerCell*nbGaussPt]);
141 INTERP_KERNEL::GaussInfo calculator(ct,gsCoods2,nbGaussPt,refCoods2,nbPtsPerCell);
142 calculator.initLocalInfo();
143 const std::vector<double>& wgths(loc.getGaussWeights());
144 for(int i=0;i<nbGaussPt;i++)
146 const double *pt0(calculator.getFunctionValues(i));
147 if(ct!=INTERP_KERNEL::NORM_HEXA27)
148 std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
151 for(int j=0;j<27;j++)
152 shape[nbPtsPerCell*i+j]=pt0[MEDMeshMultiLev::HEXA27_PERM_ARRAY[j]];
155 m[vtkType]=nbGaussPt;
156 def->Initialize(vtkType,nbPtsPerCell,nbGaussPt,shape,const_cast<double *>(&wgths[0]));
158 key->Set(elga->GetInformation(),def,vtkType);
159 key->Set(vtkd->GetInformation(),def,vtkType);
160 defs.push_back(std::pair< vtkQuadratureSchemeDefinition *, unsigned char >(def,vtkType));
163 vtkIdType ncell(ds->GetNumberOfCells());
164 int *pt(new int[ncell]),offset(0);
165 for(vtkIdType cellId=0;cellId<ncell;cellId++)
167 vtkCell *cell(ds->GetCell(cellId));
168 int delta(m[cell->GetCellType()]);
172 elga->GetInformation()->Set(MEDUtilities::ELGA(),1);
173 elga->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
174 std::ostringstream oss; oss << "ELGA" << "@" << _loc_names.size();
175 std::string ossStr(oss.str());
176 elga->SetName(ossStr.c_str());
177 elga->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
178 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elga->GetName());
179 elgas.push_back(elga);
183 _defs.push_back(defs);
187 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
189 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
190 ds->GetCellData()->AddArray(*it);
195 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
197 for(std::vector< std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > >::const_iterator it0=_defs.begin();it0!=_defs.end();it0++)
198 for(std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
199 (*it1).first->Delete();
205 class MEDFileVTKTraits
208 typedef void VtkType;
213 class MEDFileVTKTraits<int>
216 typedef vtkIntArray VtkType;
217 typedef MEDCoupling::DataArrayInt MCType;
221 class MEDFileVTKTraits<double>
224 typedef vtkDoubleArray VtkType;
225 typedef MEDCoupling::DataArrayDouble MCType;
229 void AssignDataPointerToVTK(typename MEDFileVTKTraits<T>::VtkType *vtkTab, typename MEDFileVTKTraits<T>::MCType *mcTab, bool noCpyNumNodes)
232 vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),1,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE);
234 { vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),0,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE); mcTab->accessToMemArray().setSpecificDeallocator(0); }
237 // here copy is always assumed.
238 template<class VTKT, class MCT>
239 void AssignDataPointerOther(VTKT *vtkTab, MCT *mcTab, int nbElems)
241 vtkTab->SetVoidArray(reinterpret_cast<unsigned char *>(mcTab->getPointer()),nbElems,0,VTKT::VTK_DATA_ARRAY_FREE);
242 mcTab->accessToMemArray().setSpecificDeallocator(0);
247 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
251 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>& arr):MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
253 std::vector< std::vector<MEDCoupling::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
255 throw INTERP_KERNEL::Exception("There is a big internal problem in MEDLoader ! The field time spitting has failed ! A CRASH will occur soon !");
256 if(typs[0].size()!=1)
257 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 !");
258 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDiscretization> fd(MEDCouplingFieldDiscretization::New(typs[0][0]));
259 std::ostringstream oss2; oss2 << (operator->())->getName() << ZE_SEP << fd->getRepr();
263 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
265 MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>::operator=(other);
268 _ze_name=other._ze_name;
269 _ze_full_name.clear();
273 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
278 int MEDFileFieldRepresentationLeavesArrays::getId() const
283 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
285 return _ze_full_name;
288 const char *MEDFileFieldRepresentationLeavesArrays::getZeNameC() const
290 return _ze_full_name.c_str();
293 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
295 vtkIdType refId(sil->AddChild(root,edge));
296 names.push_back(_ze_name);
298 if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
300 sil->AddChild(refId,edge);
301 names.push_back(std::string());
305 void MEDFileFieldRepresentationLeavesArrays::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
307 std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
308 _ze_full_name=oss3.str();
311 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
316 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
318 bool ret(_activated!=status);
323 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshStruct *mst, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
325 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
326 tr->setNumberOfTS((operator->())->getNumberOfTS());
328 for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
330 MCAuto<MEDFileAnyTypeField1TS> f1ts((operator->())->getTimeStepAtPos(tr->getCurrent()));
331 MEDFileAnyTypeField1TS *f1tsPtr(f1ts);
332 MEDFileField1TS *f1tsPtrDbl(dynamic_cast<MEDFileField1TS *>(f1tsPtr));
333 MEDFileIntField1TS *f1tsPtrInt(dynamic_cast<MEDFileIntField1TS *>(f1tsPtr));
334 DataArray *crudeArr(0),*postProcessedArr(0);
336 crudeArr=f1tsPtrDbl->getUndergroundDataArray();
338 crudeArr=f1tsPtrInt->getUndergroundDataArray();
340 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment !");
341 MEDFileField1TSStructItem fsst(MEDFileField1TSStructItem::BuildItemFrom(f1ts,mst));
342 f1ts->loadArraysIfNecessary();
343 MCAuto<DataArray> v(mml->buildDataArray(fsst,globs,crudeArr));
346 std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
348 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
349 vtkFieldData *att(0);
354 att=ds->GetCellData();
359 att=ds->GetPointData();
364 att=ds->GetFieldData();
369 att=ds->GetFieldData();
373 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
377 DataArray *vPtr(v); DataArrayDouble *vd(static_cast<DataArrayDouble *>(vPtr));
378 vtkDoubleArray *vtkd(vtkDoubleArray::New());
379 vtkd->SetNumberOfComponents(vd->getNumberOfComponents());
380 for(int i=0;i<vd->getNumberOfComponents();i++)
381 vtkd->SetComponentName(i,vd->getInfoOnComponent(i).c_str());
382 AssignDataPointerToVTK<double>(vtkd,vd,postProcessedArr==crudeArr);
383 std::string name(tr->buildName(f1ts->getName()));
384 vtkd->SetName(name.c_str());
387 if(discs[0]==ON_GAUSS_PT)
390 _elga_cmp.findOrCreate(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp,internalInfo);
392 if(discs[0]==ON_GAUSS_NE)
394 vtkIdTypeArray *elno(vtkIdTypeArray::New());
395 elno->SetNumberOfComponents(1);
396 vtkIdType ncell(ds->GetNumberOfCells());
397 int *pt(new int[ncell]),offset(0);
398 std::set<int> cellTypes;
399 for(vtkIdType cellId=0;cellId<ncell;cellId++)
401 vtkCell *cell(ds->GetCell(cellId));
402 int delta(cell->GetNumberOfPoints());
403 cellTypes.insert(cell->GetCellType());
407 elno->GetInformation()->Set(MEDUtilities::ELNO(),1);
408 elno->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
409 std::string nameElno("ELNO"); nameElno+="@"; nameElno+=name;
410 elno->SetName(nameElno.c_str());
411 ds->GetCellData()->AddArray(elno);
412 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elno->GetName());
413 elno->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
415 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
416 for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
418 const unsigned char *pos(std::find(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,*it));
419 if(pos==MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
421 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)std::distance(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos));
422 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
423 int nbGaussPt(cm.getNumberOfNodes()),dim(cm.getDimension());
424 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
425 double *shape(new double[nbGaussPt*nbGaussPt]);
427 const double *gsCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));//GetLocsFromGeometricType
428 const double *refCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));
429 const double *weights(MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(ct,dummy));
430 std::vector<double> gsCoords2(gsCoords,gsCoords+nbGaussPt*dim),refCoords2(refCoords,refCoords+nbGaussPt*dim);
431 INTERP_KERNEL::GaussInfo calculator(ct,gsCoords2,nbGaussPt,refCoords2,nbGaussPt);
432 calculator.initLocalInfo();
433 for(int i=0;i<nbGaussPt;i++)
435 const double *pt0(calculator.getFunctionValues(i));
436 std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
438 def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
440 key->Set(elno->GetInformation(),def,*it);
441 key->Set(vtkd->GetInformation(),def,*it);
450 DataArray *vPtr(v); DataArrayInt *vi(static_cast<DataArrayInt *>(vPtr));
451 vtkIntArray *vtkd(vtkIntArray::New());
452 vtkd->SetNumberOfComponents(vi->getNumberOfComponents());
453 for(int i=0;i<vi->getNumberOfComponents();i++)
454 vtkd->SetComponentName(i,vi->getVarOnComponent(i).c_str());
455 AssignDataPointerToVTK<int>(vtkd,vi,postProcessedArr==crudeArr);
456 std::string name(tr->buildName(f1ts->getName()));
457 vtkd->SetName(name.c_str());
462 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
466 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
468 _elga_cmp.appendELGAIfAny(ds);
473 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
477 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves(const std::vector< MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS> >& arr,
478 const MEDCoupling::MCAuto<MEDCoupling::MEDFileFastCellSupportComparator>& fsp):_arrays(arr.size()),_fsp(fsp),_cached_ds(0)
480 for(std::size_t i=0;i<arr.size();i++)
481 _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
484 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
487 _cached_ds->Delete();
490 bool MEDFileFieldRepresentationLeaves::empty() const
492 const MEDFileFastCellSupportComparator *fcscp(_fsp);
493 return fcscp==0 || _arrays.empty();
496 void MEDFileFieldRepresentationLeaves::setId(int& id) const
498 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
502 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
504 return _arrays[0]->getMeshName();
507 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
509 return (int)_arrays.size();
512 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
514 return _arrays[0]->getNumberOfTS();
517 void MEDFileFieldRepresentationLeaves::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
519 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
520 (*it).computeFullNameInLeaves(tsName,meshName,comSupStr);
524 * \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.
526 void MEDFileFieldRepresentationLeaves::feedSIL(const MEDCoupling::MEDFileMeshes *ms, const std::string& meshName, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
528 vtkIdType root2(sil->AddChild(root,edge));
529 names.push_back(std::string("Arrs"));
530 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
531 (*it).feedSIL(sil,root2,edge,names);
533 vtkIdType root3(sil->AddChild(root,edge));
534 names.push_back(std::string("InfoOnGeoType"));
535 const MEDCoupling::MEDFileMesh *m(0);
537 m=ms->getMeshWithName(meshName);
538 const MEDCoupling::MEDFileFastCellSupportComparator *fsp(_fsp);
539 if(!fsp || fsp->getNumberOfTS()==0)
541 std::vector< INTERP_KERNEL::NormalizedCellType > gts(fsp->getGeoTypesAt(0,m));
542 for(std::vector< INTERP_KERNEL::NormalizedCellType >::const_iterator it2=gts.begin();it2!=gts.end();it2++)
544 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(*it2));
545 std::string cmStr(cm.getRepr()); cmStr=cmStr.substr(5);//skip "NORM_"
546 sil->AddChild(root3,edge);
547 names.push_back(cmStr);
551 bool MEDFileFieldRepresentationLeaves::containId(int id) const
553 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
554 if((*it).getId()==id)
559 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
561 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
562 if((*it).getZeName()==name)
570 void MEDFileFieldRepresentationLeaves::dumpState(std::map<std::string,bool>& status) const
572 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
573 status[(*it).getZeName()]=(*it).getStatus();
576 bool MEDFileFieldRepresentationLeaves::isActivated() const
578 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
579 if((*it).getStatus())
584 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
586 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
588 os << " - " << (*it0).getZeName() << " (";
589 if((*it0).getStatus())
593 os << ")" << std::endl;
597 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
599 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
600 (*it).setStatus(true);
603 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
605 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
606 if((*it).getId()==id)
608 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
611 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
614 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
615 std::vector<double> ret;
616 std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
617 return tk.getTimeStepsRegardingPolicy(dtits,ret);
620 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
623 return _arrays[0]->getTimeSteps(ts);
627 return std::vector< std::pair<int,int> >();
631 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
633 std::ostringstream oss;
634 oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
635 std::vector<double> ret1;
636 std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
637 std::size_t sz(ret1.size());
638 for(std::size_t i=0;i<sz;i++)
640 oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
643 std::string tmp(oss.str());
644 if(tmp.size()>200 && i!=sz-1)
654 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshes *meshes, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
657 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
658 MCAuto<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
659 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
660 if((*it).getStatus())
662 (*it).appendFields(tr,globs,mml,mst,ds,internalInfo);
663 (*it).appendELGAIfAny(ds);
667 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
669 DataArrayDouble *coordsMC(0);
670 DataArrayByte *typesMC(0);
671 DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
672 bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
673 MCAuto<DataArrayDouble> coordsSafe(coordsMC);
674 MCAuto<DataArrayByte> typesSafe(typesMC);
675 MCAuto<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
677 int nbOfCells(typesSafe->getNbOfElems());
678 vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
679 vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
680 AssignDataPointerOther<vtkUnsignedCharArray,DataArrayByte>(cellTypes,typesSafe,nbOfCells);
681 vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
682 AssignDataPointerOther<vtkIdTypeArray,DataArrayInt>(cellLocations,cellLocationsSafe,nbOfCells);
683 vtkCellArray *cells(vtkCellArray::New());
684 vtkIdTypeArray *cells2(vtkIdTypeArray::New());
685 AssignDataPointerOther<vtkIdTypeArray,DataArrayInt>(cells2,cellsSafe,cellsSafe->getNbOfElems());
686 cells->SetCells(nbOfCells,cells2);
688 if(faceLocationsMC!=0 && facesMC!=0)
690 vtkIdTypeArray *faces(vtkIdTypeArray::New());
691 AssignDataPointerOther<vtkIdTypeArray,DataArrayInt>(faces,facesSafe,facesSafe->getNbOfElems());
692 vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
693 AssignDataPointerOther<vtkIdTypeArray,DataArrayInt>(faceLocations,faceLocationsSafe,faceLocationsSafe->getNbOfElems());
694 ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
695 faceLocations->Delete();
699 ret->SetCells(cellTypes,cellLocations,cells);
701 cellLocations->Delete();
703 vtkPoints *pts(vtkPoints::New());
704 vtkDoubleArray *pts2(vtkDoubleArray::New());
705 pts2->SetNumberOfComponents(3);
706 AssignDataPointerToVTK<double>(pts2,coordsSafe,statusOfCoords);
715 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(MEDCoupling::MEDCMeshMultiLev *mm) const
718 std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
719 vtkDoubleArray *vtkTmp(0);
720 vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
721 std::size_t dim(arrs.size());
723 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
724 int sizePerAxe[3]={1,1,1};
725 sizePerAxe[0]=arrs[0]->getNbOfElems();
727 sizePerAxe[1]=arrs[1]->getNbOfElems();
729 sizePerAxe[2]=arrs[2]->getNbOfElems();
730 ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
731 vtkTmp=vtkDoubleArray::New();
732 vtkTmp->SetNumberOfComponents(1);
733 AssignDataPointerToVTK<double>(vtkTmp,arrs[0],isInternal);
734 ret->SetXCoordinates(vtkTmp);
739 vtkTmp=vtkDoubleArray::New();
740 vtkTmp->SetNumberOfComponents(1);
741 AssignDataPointerToVTK<double>(vtkTmp,arrs[1],isInternal);
742 ret->SetYCoordinates(vtkTmp);
748 vtkTmp=vtkDoubleArray::New();
749 vtkTmp->SetNumberOfComponents(1);
750 AssignDataPointerToVTK<double>(vtkTmp,arrs[2],isInternal);
751 ret->SetZCoordinates(vtkTmp);
758 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(MEDCoupling::MEDCurveLinearMeshMultiLev *mm) const
760 int meshStr[3]={1,1,1};
761 DataArrayDouble *coords(0);
762 std::vector<int> nodeStrct;
764 mm->buildVTUArrays(coords,nodeStrct,isInternal);
765 std::size_t dim(nodeStrct.size());
767 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
768 meshStr[0]=nodeStrct[0];
770 meshStr[1]=nodeStrct[1];
772 meshStr[2]=nodeStrct[2];
773 vtkStructuredGrid *ret(vtkStructuredGrid::New());
774 ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
775 vtkDoubleArray *da(vtkDoubleArray::New());
776 da->SetNumberOfComponents(3);
777 if(coords->getNumberOfComponents()==3)
778 AssignDataPointerToVTK<double>(da,coords,isInternal);//if isIntenal==True VTK has not the ownership of double * because MEDLoader main struct has it !
781 MCAuto<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
782 AssignDataPointerToVTK<double>(da,coords2,false);//let VTK deal with double *
785 vtkPoints *points=vtkPoints::New();
786 ret->SetPoints(points);
793 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, ExportedTinyInfo *internalInfo) const
796 //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
797 MCAuto<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change
798 MCAuto<MEDMeshMultiLev> mml2(mml->prepare());
799 MEDMeshMultiLev *ptMML2(mml2);
802 MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
803 MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
804 MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
808 ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
812 ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
816 ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
819 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
820 _cached_ds=ret->NewInstance();
821 _cached_ds->ShallowCopy(ret);
825 ret=_cached_ds->NewInstance();
826 ret->ShallowCopy(_cached_ds);
829 appendFields(tr,globs,mml,meshes,ret,internalInfo);
830 // The arrays links to mesh
831 DataArrayInt *famCells(0),*numCells(0);
832 bool noCpyFamCells(false),noCpyNumCells(false);
833 ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
836 vtkIntArray *vtkTab(vtkIntArray::New());
837 vtkTab->SetNumberOfComponents(1);
838 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
839 AssignDataPointerToVTK<int>(vtkTab,famCells,noCpyFamCells);
840 ret->GetCellData()->AddArray(vtkTab);
844 ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
847 vtkIntArray *vtkTab(vtkIntArray::New());
848 vtkTab->SetNumberOfComponents(1);
849 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
850 AssignDataPointerToVTK<int>(vtkTab,numCells,noCpyNumCells);
851 ret->GetCellData()->AddArray(vtkTab);
855 // The arrays links to mesh
856 DataArrayInt *famNodes(0),*numNodes(0);
857 bool noCpyFamNodes(false),noCpyNumNodes(false);
858 ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
861 vtkIntArray *vtkTab(vtkIntArray::New());
862 vtkTab->SetNumberOfComponents(1);
863 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
864 AssignDataPointerToVTK<int>(vtkTab,famNodes,noCpyFamNodes);
865 ret->GetPointData()->AddArray(vtkTab);
869 ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
872 vtkIntArray *vtkTab(vtkIntArray::New());
873 vtkTab->SetNumberOfComponents(1);
874 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
875 AssignDataPointerToVTK<int>(vtkTab,numNodes,noCpyNumNodes);
876 ret->GetPointData()->AddArray(vtkTab);
880 // Global Node Ids if any ! (In // mode)
881 DataArrayInt *gni(ptMML2->retrieveGlobalNodeIdsIfAny());
884 vtkIntArray *vtkTab(vtkIntArray::New());
885 vtkTab->SetNumberOfComponents(1);
886 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::GLOBAL_NODE_ID_NAME);
887 AssignDataPointerToVTK<int>(vtkTab,gni,false);
888 ret->GetPointData()->AddArray(vtkTab);
895 //////////////////////
897 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
901 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
904 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
905 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
906 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
907 ret+=(*it2).getNumberOfArrays();
911 void MEDFileFieldRepresentationTree::assignIds() const
914 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
915 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
916 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
920 void MEDFileFieldRepresentationTree::computeFullNameInLeaves() const
922 std::size_t it0Cnt(0);
923 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
925 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
926 std::string tsName(oss.str());
927 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
929 std::string meshName((*it1)[0].getMeshName());
930 std::size_t it2Cnt(0);
931 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
933 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
934 std::string comSupStr(oss2.str());
935 (*it2).computeFullNameInLeaves(tsName,meshName,comSupStr);
941 void MEDFileFieldRepresentationTree::activateTheFirst() const
943 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
944 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
945 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
947 (*it2).activateAllArrays();
952 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
954 std::size_t it0Cnt(0);
955 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
957 vtkIdType InfoOnTSId(sil->AddChild(root,edge));
958 names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
960 vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
961 std::vector<double> ts;
962 std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
963 std::size_t nbOfTS(dtits.size());
964 std::ostringstream oss3; oss3 << nbOfTS;
965 names.push_back(oss3.str());
966 for(std::size_t i=0;i<nbOfTS;i++)
968 std::ostringstream oss4; oss4 << dtits[i].first;
969 vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
970 names.push_back(oss4.str());
971 std::ostringstream oss5; oss5 << dtits[i].second;
972 vtkIdType ItId(sil->AddChild(DtId,edge));
973 names.push_back(oss5.str());
974 std::ostringstream oss6; oss6 << ts[i];
975 sil->AddChild(ItId,edge);
976 names.push_back(oss6.str());
979 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
980 std::string tsName(oss.str());
981 vtkIdType typeId0(sil->AddChild(root,edge));
982 names.push_back(tsName);
983 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
985 std::string meshName((*it1)[0].getMeshName());
986 vtkIdType typeId1(sil->AddChild(typeId0,edge));
987 names.push_back(meshName);
988 std::size_t it2Cnt(0);
989 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
991 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
992 std::string comSupStr(oss2.str());
993 vtkIdType typeId2(sil->AddChild(typeId1,edge));
994 names.push_back(comSupStr);
995 (*it2).feedSIL(_ms,meshName,sil,typeId2,edge,names);
1001 std::string MEDFileFieldRepresentationTree::getActiveMeshName() const
1003 int dummy0(0),dummy1(0),dummy2(0);
1004 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
1005 return leaf.getMeshName();
1008 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
1010 int dummy0(0),dummy1(0),dummy2(0);
1011 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
1012 std::string ret(leaf.getMeshName());
1015 for(;i<_ms->getNumberOfMeshes();i++)
1017 m=_ms->getMeshAtPos(i);
1018 if(m->getName()==ret)
1021 if(i==_ms->getNumberOfMeshes())
1022 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
1023 vtkIdType typeId0(sil->AddChild(root,edge));
1024 names.push_back(m->getName());
1026 vtkIdType typeId1(sil->AddChild(typeId0,edge));
1027 names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
1028 std::vector<std::string> grps(m->getGroupsNames());
1029 for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
1031 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1032 names.push_back(*it0);
1033 std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
1034 for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
1036 sil->AddChild(typeId2,edge);
1037 names.push_back((*it1).c_str());
1041 vtkIdType typeId11(sil->AddChild(typeId0,edge));
1042 names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
1043 std::vector<std::string> fams(m->getFamiliesNames());
1044 for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
1046 sil->AddChild(typeId11,edge);
1047 int famId(m->getFamilyId((*it00).c_str()));
1048 std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
1049 names.push_back(oss.str());
1054 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
1056 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1057 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1058 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1059 if((*it2).containId(id))
1060 return (*it2).getLeafArr(id);
1061 throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
1064 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
1066 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1067 return elt.getZeName();
1070 const char *MEDFileFieldRepresentationTree::getNameOfC(int id) const
1072 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1073 return elt.getZeNameC();
1076 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1078 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1079 return elt.getStatus();
1082 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1085 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1086 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1087 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1088 if((*it2).containZeName(name,ret))
1090 std::ostringstream msg; msg << "MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name \"" << name << "\" !";
1091 throw INTERP_KERNEL::Exception(msg.str().c_str());
1094 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1096 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1097 bool ret(elt.setStatus(status));//to be implemented
1101 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1104 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1105 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1106 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1107 ret=std::max(ret,(*it2).getNumberOfTS());
1114 void MEDFileFieldRepresentationTree::loadInMemory(MEDCoupling::MEDFileFields *fields, MEDCoupling::MEDFileMeshes *meshes)
1116 _fields=fields; _ms=meshes;
1117 if(_fields.isNotNull())
1122 if(_fields.isNull())
1124 _fields=BuildFieldFromMeshes(_ms);
1128 AppendFieldFromMeshes(_ms,_fields);
1130 _ms->cartesianizeMe();
1131 _fields->removeFieldsWithoutAnyTimeStep();
1132 std::vector<std::string> meshNames(_ms->getMeshesNames());
1133 std::vector< MCAuto<MEDFileFields> > fields_per_mesh(meshNames.size());
1134 for(std::size_t i=0;i<meshNames.size();i++)
1136 fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1138 std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1140 for(std::vector< MCAuto<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1142 for(int j=0;j<(*fields)->getNumberOfFields();j++)
1144 MCAuto<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1145 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1147 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::const_iterator it=tmp.begin();it!=tmp.end();it++)
1149 if(!(*it)->presenceOfMultiDiscPerGeoType())
1150 allFMTSLeavesToDisplaySafe.push_back(*it);
1152 {// The case of some parts of field have more than one discretization per geo type.
1153 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > subTmp((*it)->splitMultiDiscrPerGeoTypes());
1154 std::size_t it0Cnt(0);
1155 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::iterator it0=subTmp.begin();it0!=subTmp.end();it0++,it0Cnt++)//not const because setName
1157 std::ostringstream oss; oss << (*it0)->getName() << "_" << std::setfill('M') << std::setw(3) << it0Cnt;
1158 (*it0)->setName(oss.str());
1159 allFMTSLeavesToDisplaySafe.push_back(*it0);
1166 std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1167 for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1169 allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1171 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1172 // memory safety part
1173 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1174 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1176 allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1177 for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1179 allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1180 allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1183 // end of memory safety part
1184 // 1st : timesteps, 2nd : meshName, 3rd : common support
1185 this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1186 for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1188 std::vector< std::string > meshNamesLoc;
1189 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1190 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1192 std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1193 std::vector< std::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1194 if(it==meshNamesLoc.end())
1196 meshNamesLoc.push_back(meshName);
1197 splitByMeshName.resize(splitByMeshName.size()+1);
1198 splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1201 splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1203 _data_structure[i].resize(meshNamesLoc.size());
1204 for(std::size_t j=0;j<splitByMeshName.size();j++)
1206 std::vector< MCAuto<MEDFileFastCellSupportComparator> > fsp;
1207 std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1208 for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1209 sbmn[k]=splitByMeshName[j][k];
1210 //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1211 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1212 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1213 this->_data_structure[i][j].resize(commonSupSplit.size());
1214 for(std::size_t k=0;k<commonSupSplit.size();k++)
1216 commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1217 for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1219 commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1220 commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1223 for(std::size_t k=0;k<commonSupSplit.size();k++)
1224 this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1227 this->removeEmptyLeaves();
1229 this->computeFullNameInLeaves();
1232 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts)
1234 MCAuto<MEDFileMeshes> ms;
1235 MCAuto<MEDFileFields> fields;
1238 if((iPart==-1 && nbOfParts==-1) || (iPart==0 && nbOfParts==1))
1240 MCAuto<MEDFileMeshSupports> msups(MEDFileMeshSupports::New(fileName));
1241 MCAuto<MEDFileStructureElements> mse(MEDFileStructureElements::New(fileName,msups));
1242 ms=MEDFileMeshes::New(fileName);
1243 fields=MEDFileFields::NewWithDynGT(fileName,mse,false);//false is important to not read the values
1244 if(ms->presenceOfStructureElements())
1246 fields->loadArrays();
1247 MEDFileBlowStrEltUp::DealWithSE(fields,ms,mse);
1249 int nbMeshes(ms->getNumberOfMeshes());
1250 for(int i=0;i<nbMeshes;i++)
1252 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1253 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1255 tmp2->forceComputationOfParts();
1260 #ifdef MEDREADER_USE_MPI
1261 ms=ParaMEDFileMeshes::New(iPart,nbOfParts,fileName);
1262 int nbMeshes(ms->getNumberOfMeshes());
1263 for(int i=0;i<nbMeshes;i++)
1265 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1266 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1268 MCAuto<DataArrayInt> tmp3(tmp2->zipCoords());
1270 fields=MEDFileFields::LoadPartOf(fileName,false,ms);//false is important to not read the values
1272 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::loadMainStructureOfFile : request for iPart/nbOfParts=" << iPart << "/" << nbOfParts << " whereas Plugin not compiled with MPI !";
1273 throw INTERP_KERNEL::Exception(oss.str().c_str());
1279 MCAuto<MEDCoupling::SauvReader> sr(MEDCoupling::SauvReader::New(fileName));
1280 MCAuto<MEDCoupling::MEDFileData> mfd(sr->loadInMEDFileDS());
1281 ms=mfd->getMeshes(); ms->incrRef();
1282 int nbMeshes(ms->getNumberOfMeshes());
1283 for(int i=0;i<nbMeshes;i++)
1285 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1286 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1288 tmp2->forceComputationOfParts();
1290 fields=mfd->getFields();
1291 if(fields.isNotNull())
1294 loadInMemory(fields,ms);
1297 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1299 std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1300 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1302 std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1303 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1305 std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1306 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1308 newSD1.push_back(*it2);
1310 newSD0.push_back(newSD1);
1313 newSD.push_back(newSD0);
1317 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1319 if(compInfos.size()!=1)
1321 return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1324 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1326 return _data_structure[0][0][0].getMeshName();
1329 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1332 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1333 return leaf.getTimeSteps(tk);
1336 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, ExportedTinyInfo *internalInfo) const
1339 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1340 meshName=leaf.getMeshName();
1341 std::vector<double> ts(leaf.getTimeSteps(tk));
1342 std::size_t zeTimeId(0);
1345 std::vector<double> ts2(ts.size());
1346 std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1347 std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1348 zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1351 if(zeTimeId==(int)ts.size())
1352 zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1353 if(zeTimeId==(int)ts.size())
1354 {//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.
1355 //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1357 double valAttachedToPos(-std::numeric_limits<double>::max());
1358 for(std::size_t i=0;i<ts.size();i++)
1362 if(ts[i]>valAttachedToPos)
1365 valAttachedToPos=ts[i];
1370 {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1371 valAttachedToPos=std::numeric_limits<double>::max();
1372 for(std::size_t i=0;i<ts.size();i++)
1374 if(ts[i]<valAttachedToPos)
1377 valAttachedToPos=ts[i];
1382 std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1383 std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1384 oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1385 std::cerr << oss.str() << std::endl;
1389 tr=new MEDStdTimeReq((int)zeTimeId);
1391 tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
1392 vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms,internalInfo));
1397 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1399 int nbOfActivated(0);
1400 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1401 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1402 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1403 if((*it2).isActivated())
1405 if(nbOfActivated!=1)
1407 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1408 throw INTERP_KERNEL::Exception(oss.str().c_str());
1410 int i0(0),i1(0),i2(0);
1411 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1412 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1413 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1414 if((*it2).isActivated())
1416 lev0=i0; lev1=i1; lev2=i2;
1419 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1422 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1425 os << "#############################################" << std::endl;
1426 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1429 os << "TS" << i << std::endl;
1430 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1433 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1436 os << " " << (*it2).getMeshName() << std::endl;
1437 os << " Comp" << k << std::endl;
1438 (*it2).printMySelf(os);
1442 os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1445 std::map<std::string,bool> MEDFileFieldRepresentationTree::dumpState() const
1447 std::map<std::string,bool> ret;
1448 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1449 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1450 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1451 (*it2).dumpState(ret);
1455 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms, MEDCoupling::MEDFileFields *ret)
1458 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::AppendFieldFromMeshes : internal error ! NULL ret !");
1459 for(int i=0;i<ms->getNumberOfMeshes();i++)
1461 MEDFileMesh *mm(ms->getMeshAtPos(i));
1462 std::vector<int> levs(mm->getNonEmptyLevels());
1463 MEDCoupling::MCAuto<MEDCoupling::MEDFileField1TS> f1tsMultiLev(MEDCoupling::MEDFileField1TS::New());
1464 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1467 for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1469 std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1470 for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1472 MEDCoupling::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1473 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1475 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1476 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1479 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1480 f1tsMultiLev->setFieldNoProfileSBT(f);
1485 std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1486 if(levsExt.size()==levs.size()+1)
1488 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(1));
1489 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_NODES));
1491 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1492 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1493 arr->iota(); f->setArray(arr);
1494 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1495 f1tsMultiLev->setFieldNoProfileSBT(f);
1503 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(0));
1504 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1506 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1507 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1510 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1511 f1tsMultiLev->setFieldNoProfileSBT(f);
1514 MEDCoupling::MCAuto<MEDCoupling::MEDFileFieldMultiTS> fmtsMultiLev(MEDCoupling::MEDFileFieldMultiTS::New());
1515 fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1516 ret->pushField(fmtsMultiLev);
1520 std::string MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh(const std::string& meshName, const MEDCoupling::MEDFileFields *ret)
1522 const char KEY_STR_TO_AVOID_COLLIDE[]="MESH@";
1524 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh : internal error ! NULL ret !");
1525 std::vector<std::string> fieldNamesAlreadyExisting(ret->getFieldsNames());
1526 if(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),meshName)==fieldNamesAlreadyExisting.end())
1528 std::string tmpName(KEY_STR_TO_AVOID_COLLIDE); tmpName+=meshName;
1529 while(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),tmpName)!=fieldNamesAlreadyExisting.end())
1530 tmpName=std::string(KEY_STR_TO_AVOID_COLLIDE)+tmpName;
1534 MEDCoupling::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms)
1536 MEDCoupling::MCAuto<MEDCoupling::MEDFileFields> ret(MEDCoupling::MEDFileFields::New());
1537 AppendFieldFromMeshes(ms,ret);
1541 std::vector<std::string> MEDFileFieldRepresentationTree::SplitFieldNameIntoParts(const std::string& fullFieldName, char sep)
1543 std::vector<std::string> ret;
1545 while(pos!=std::string::npos)
1547 std::size_t curPos(fullFieldName.find_first_of(sep,pos));
1548 std::string elt(fullFieldName.substr(pos,curPos!=std::string::npos?curPos-pos:std::string::npos));
1550 pos=fullFieldName.find_first_not_of(sep,curPos);
1556 * Here the non regression tests.
1557 * const char inp0[]="";
1558 * const char exp0[]="";
1559 * const char inp1[]="field";
1560 * const char exp1[]="field";
1561 * const char inp2[]="_________";
1562 * const char exp2[]="_________";
1563 * const char inp3[]="field_p";
1564 * const char exp3[]="field_p";
1565 * const char inp4[]="field__p";
1566 * const char exp4[]="field_p";
1567 * const char inp5[]="field_p__";
1568 * const char exp5[]="field_p";
1569 * const char inp6[]="field_p_";
1570 * const char exp6[]="field_p";
1571 * const char inp7[]="field_____EDFGEG//sdkjf_____PP_______________";
1572 * const char exp7[]="field_EDFGEG//sdkjf_PP";
1573 * const char inp8[]="field_____EDFGEG//sdkjf_____PP";
1574 * const char exp8[]="field_EDFGEG//sdkjf_PP";
1575 * const char inp9[]="_field_____EDFGEG//sdkjf_____PP_______________";
1576 * const char exp9[]="field_EDFGEG//sdkjf_PP";
1577 * const char inp10[]="___field_____EDFGEG//sdkjf_____PP_______________";
1578 * const char exp10[]="field_EDFGEG//sdkjf_PP";
1580 std::string MEDFileFieldRepresentationTree::PostProcessFieldName(const std::string& fullFieldName)
1582 const char SEP('_');
1583 std::vector<std::string> v(SplitFieldNameIntoParts(fullFieldName,SEP));
1585 return fullFieldName;//should never happen
1589 return fullFieldName;
1593 std::string ret(v[0]);
1594 for(std::size_t i=1;i<v.size();i++)
1599 { ret+=SEP; ret+=v[i]; }
1605 return fullFieldName;
1611 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1615 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1620 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1622 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1624 throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1630 * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1631 * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1632 * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1634 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1636 std::size_t sz(ts.size());
1637 bool isInHumanRange(true);
1639 for(std::size_t i=0;i<sz;i++)
1642 if(ts[i]<=-1e299 || ts[i]>=1e299)
1643 isInHumanRange=false;
1646 return processedUsingPairOfIds(tsPairs);
1648 return processedUsingPairOfIds(tsPairs);
1649 _postprocessed_time=ts;
1650 return getPostProcessedTime();
1655 * idem than 0, except that ts is preaccumulated before invoking policy 0.
1657 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1659 std::size_t sz(ts.size());
1660 std::vector<double> ts2(sz);
1662 for(std::size_t i=0;i<sz;i++)
1667 return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1670 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1672 std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1676 void TimeKeeper::printSelf(std::ostream& oss) const
1678 std::size_t sz(_activated_ts.size());
1679 for(std::size_t i=0;i<sz;i++)
1681 oss << "(" << i << "," << _activated_ts[i].first << "), ";
1685 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1687 std::size_t sz(_activated_ts.size());
1688 std::vector<bool> ret(sz);
1689 for(std::size_t i=0;i<sz;i++)
1691 ret[i]=_activated_ts[i].first;
1696 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1698 std::size_t sz(tsPairs.size());
1699 std::set<int> s0,s1;
1700 for(std::size_t i=0;i<sz;i++)
1701 { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1704 _postprocessed_time.resize(sz);
1705 for(std::size_t i=0;i<sz;i++)
1706 _postprocessed_time[i]=(double)tsPairs[i].first;
1707 return getPostProcessedTime();
1711 _postprocessed_time.resize(sz);
1712 for(std::size_t i=0;i<sz;i++)
1713 _postprocessed_time[i]=(double)tsPairs[i].second;
1714 return getPostProcessedTime();
1716 //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1717 _postprocessed_time.resize(sz);
1718 for(std::size_t i=0;i<sz;i++)
1719 _postprocessed_time[i]=(double)i;
1720 return getPostProcessedTime();
1723 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1725 _activated_ts.resize(maxNumberOfTS);
1726 for(int i=0;i<maxNumberOfTS;i++)
1728 std::ostringstream oss; oss << "000" << i;
1729 _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());