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 "vtkMutableDirectedGraph.h"
54 #include "vtkDataArrayTemplate.h"
56 using namespace MEDCoupling;
58 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
60 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
62 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
64 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
66 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
68 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
70 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
72 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
74 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
76 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
78 vtkIdTypeArray *ELGACmp::findOrCreate(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds, bool& isNew) const
80 vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
89 return createNew(globs,locsReallyUsed,vtkd,ds);
93 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd) const
95 std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
96 if(it==_loc_names.end())
98 std::size_t pos(std::distance(_loc_names.begin(),it));
99 vtkIdTypeArray *ret(_elgas[pos]);
100 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
101 for(std::vector<std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it=_defs[pos].begin();it!=_defs[pos].end();it++)
103 key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
105 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
109 vtkIdTypeArray *ELGACmp::createNew(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDoubleArray *vtkd, vtkDataSet *ds) const
111 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
112 std::vector< std::vector<std::string> > locNames(_loc_names);
113 std::vector<vtkIdTypeArray *> elgas(_elgas);
114 std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > defs;
116 std::vector< std::vector<std::string> >::const_iterator it(std::find(locNames.begin(),locNames.end(),locsReallyUsed));
117 if(it!=locNames.end())
118 throw INTERP_KERNEL::Exception("ELGACmp::createNew : Method is expected to be called after isExisting call ! Entry already exists !");
119 locNames.push_back(locsReallyUsed);
120 vtkIdTypeArray *elga(vtkIdTypeArray::New());
121 elga->SetNumberOfComponents(1);
122 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
123 std::map<unsigned char,int> m;
124 for(std::vector<std::string>::const_iterator it=locsReallyUsed.begin();it!=locsReallyUsed.end();it++)
126 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
127 const MEDFileFieldLoc& loc(globs->getLocalization((*it).c_str()));
128 INTERP_KERNEL::NormalizedCellType ct(loc.getGeoType());
129 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
130 int nbGaussPt(loc.getNbOfGaussPtPerCell()),nbPtsPerCell((int)cm.getNumberOfNodes()),dimLoc(loc.getDimension());
131 // 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.
132 std::vector<double> gsCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getGaussCoords()));
133 std::vector<double> refCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getRefCoords()));
134 double *shape(new double[nbPtsPerCell*nbGaussPt]);
135 INTERP_KERNEL::GaussInfo calculator(ct,gsCoods2,nbGaussPt,refCoods2,nbPtsPerCell);
136 calculator.initLocalInfo();
137 const std::vector<double>& wgths(loc.getGaussWeights());
138 for(int i=0;i<nbGaussPt;i++)
140 const double *pt0(calculator.getFunctionValues(i));
141 if(ct!=INTERP_KERNEL::NORM_HEXA27)
142 std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
145 for(int j=0;j<27;j++)
146 shape[nbPtsPerCell*i+j]=pt0[MEDMeshMultiLev::HEXA27_PERM_ARRAY[j]];
149 unsigned char vtkType(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE[ct]);
150 m[vtkType]=nbGaussPt;
151 def->Initialize(vtkType,nbPtsPerCell,nbGaussPt,shape,const_cast<double *>(&wgths[0]));
153 key->Set(elga->GetInformation(),def,vtkType);
154 key->Set(vtkd->GetInformation(),def,vtkType);
155 defs.push_back(std::pair< vtkQuadratureSchemeDefinition *, unsigned char >(def,vtkType));
158 vtkIdType ncell(ds->GetNumberOfCells());
159 int *pt(new int[ncell]),offset(0);
160 for(vtkIdType cellId=0;cellId<ncell;cellId++)
162 vtkCell *cell(ds->GetCell(cellId));
163 int delta(m[cell->GetCellType()]);
167 elga->GetInformation()->Set(MEDUtilities::ELGA(),1);
168 elga->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
169 std::ostringstream oss; oss << "ELGA" << "@" << _loc_names.size();
170 std::string ossStr(oss.str());
171 elga->SetName(ossStr.c_str());
172 elga->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
173 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elga->GetName());
174 elgas.push_back(elga);
178 _defs.push_back(defs);
181 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
183 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
184 ds->GetCellData()->AddArray(*it);
189 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
191 for(std::vector< std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > >::const_iterator it0=_defs.begin();it0!=_defs.end();it0++)
192 for(std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
193 (*it1).first->Delete();
198 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
202 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>& arr):MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
204 std::vector< std::vector<MEDCoupling::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
206 throw INTERP_KERNEL::Exception("There is a big internal problem in MEDLoader ! The field time spitting has failed ! A CRASH will occur soon !");
207 if(typs[0].size()!=1)
208 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 !");
209 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDiscretization> fd(MEDCouplingFieldDiscretization::New(typs[0][0]));
210 std::ostringstream oss2; oss2 << (operator->())->getName() << ZE_SEP << fd->getRepr();
214 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
216 MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>::operator=(other);
219 _ze_name=other._ze_name;
220 _ze_full_name.clear();
224 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
229 int MEDFileFieldRepresentationLeavesArrays::getId() const
234 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
236 return _ze_full_name;
239 const char *MEDFileFieldRepresentationLeavesArrays::getZeNameC() const
241 return _ze_full_name.c_str();
244 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
246 vtkIdType refId(sil->AddChild(root,edge));
247 names.push_back(_ze_name);
249 if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
251 sil->AddChild(refId,edge);
252 names.push_back(std::string());
256 void MEDFileFieldRepresentationLeavesArrays::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
258 std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
259 _ze_full_name=oss3.str();
262 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
267 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
269 bool ret(_activated!=status);
274 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshStruct *mst, vtkDataSet *ds) const
276 const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
277 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
278 tr->setNumberOfTS((operator->())->getNumberOfTS());
280 for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
282 MCAuto<MEDFileAnyTypeField1TS> f1ts((operator->())->getTimeStepAtPos(tr->getCurrent()));
283 MEDFileAnyTypeField1TS *f1tsPtr(f1ts);
284 MEDFileField1TS *f1tsPtrDbl(dynamic_cast<MEDFileField1TS *>(f1tsPtr));
285 MEDFileIntField1TS *f1tsPtrInt(dynamic_cast<MEDFileIntField1TS *>(f1tsPtr));
286 DataArray *crudeArr(0),*postProcessedArr(0);
288 crudeArr=f1tsPtrDbl->getUndergroundDataArray();
290 crudeArr=f1tsPtrInt->getUndergroundDataArray();
292 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment !");
293 MEDFileField1TSStructItem fsst(MEDFileField1TSStructItem::BuildItemFrom(f1ts,mst));
294 f1ts->loadArraysIfNecessary();
295 MCAuto<DataArray> v(mml->buildDataArray(fsst,globs,crudeArr));
298 std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
300 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
301 vtkFieldData *att(0);
306 att=ds->GetCellData();
311 att=ds->GetPointData();
316 att=ds->GetFieldData();
321 att=ds->GetFieldData();
325 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
329 DataArray *vPtr(v); DataArrayDouble *vd(static_cast<DataArrayDouble *>(vPtr));
330 vtkDoubleArray *vtkd(vtkDoubleArray::New());
331 vtkd->SetNumberOfComponents(vd->getNumberOfComponents());
332 for(int i=0;i<vd->getNumberOfComponents();i++)
333 vtkd->SetComponentName(i,vd->getInfoOnComponent(i).c_str());
334 if(postProcessedArr!=crudeArr)
336 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vd->accessToMemArray().setSpecificDeallocator(0);
340 vtkd->SetArray(vd->getPointer(),vd->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
342 std::string name(tr->buildName(f1ts->getName()));
343 vtkd->SetName(name.c_str());
346 if(discs[0]==ON_GAUSS_PT)
349 _elga_cmp.findOrCreate(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp);
351 if(discs[0]==ON_GAUSS_NE)
353 vtkIdTypeArray *elno(vtkIdTypeArray::New());
354 elno->SetNumberOfComponents(1);
355 vtkIdType ncell(ds->GetNumberOfCells());
356 int *pt(new int[ncell]),offset(0);
357 std::set<int> cellTypes;
358 for(vtkIdType cellId=0;cellId<ncell;cellId++)
360 vtkCell *cell(ds->GetCell(cellId));
361 int delta(cell->GetNumberOfPoints());
362 cellTypes.insert(cell->GetCellType());
366 elno->GetInformation()->Set(MEDUtilities::ELNO(),1);
367 elno->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
368 std::string nameElno("ELNO"); nameElno+="@"; nameElno+=name;
369 elno->SetName(nameElno.c_str());
370 ds->GetCellData()->AddArray(elno);
371 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elno->GetName());
372 elno->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
374 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
375 for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
377 const unsigned char *pos(std::find(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,*it));
378 if(pos==MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
380 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)std::distance(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos));
381 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
382 int nbGaussPt(cm.getNumberOfNodes()),dim(cm.getDimension());
383 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
384 double *shape(new double[nbGaussPt*nbGaussPt]);
386 const double *gsCoords(MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType(ct,dummy));
387 const double *refCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));
388 const double *weights(MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(ct,dummy));
389 std::vector<double> gsCoords2(gsCoords,gsCoords+nbGaussPt*dim),refCoords2(refCoords,refCoords+nbGaussPt*dim);
390 INTERP_KERNEL::GaussInfo calculator(ct,gsCoords2,nbGaussPt,refCoords2,nbGaussPt);
391 calculator.initLocalInfo();
392 for(int i=0;i<nbGaussPt;i++)
394 const double *pt0(calculator.getFunctionValues(i));
395 std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
397 def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
399 key->Set(elno->GetInformation(),def,*it);
400 key->Set(vtkd->GetInformation(),def,*it);
409 DataArray *vPtr(v); DataArrayInt *vi(static_cast<DataArrayInt *>(vPtr));
410 vtkIntArray *vtkd(vtkIntArray::New());
411 vtkd->SetNumberOfComponents(vi->getNumberOfComponents());
412 for(int i=0;i<vi->getNumberOfComponents();i++)
413 vtkd->SetComponentName(i,vi->getVarOnComponent(i).c_str());
414 if(postProcessedArr!=crudeArr)
416 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); vi->accessToMemArray().setSpecificDeallocator(0);
420 vtkd->SetArray(vi->getPointer(),vi->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
422 std::string name(tr->buildName(f1ts->getName()));
423 vtkd->SetName(name.c_str());
428 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
432 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
434 _elga_cmp.appendELGAIfAny(ds);
439 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
443 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves(const std::vector< MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS> >& arr,
444 const MEDCoupling::MCAuto<MEDCoupling::MEDFileFastCellSupportComparator>& fsp):_arrays(arr.size()),_fsp(fsp),_cached_ds(0)
446 for(std::size_t i=0;i<arr.size();i++)
447 _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
450 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
453 _cached_ds->Delete();
456 bool MEDFileFieldRepresentationLeaves::empty() const
458 const MEDFileFastCellSupportComparator *fcscp(_fsp);
459 return fcscp==0 || _arrays.empty();
462 void MEDFileFieldRepresentationLeaves::setId(int& id) const
464 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
468 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
470 return _arrays[0]->getMeshName();
473 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
475 return (int)_arrays.size();
478 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
480 return _arrays[0]->getNumberOfTS();
483 void MEDFileFieldRepresentationLeaves::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
485 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
486 (*it).computeFullNameInLeaves(tsName,meshName,comSupStr);
490 * \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.
492 void MEDFileFieldRepresentationLeaves::feedSIL(const MEDCoupling::MEDFileMeshes *ms, const std::string& meshName, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
494 vtkIdType root2(sil->AddChild(root,edge));
495 names.push_back(std::string("Arrs"));
496 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
497 (*it).feedSIL(sil,root2,edge,names);
499 vtkIdType root3(sil->AddChild(root,edge));
500 names.push_back(std::string("InfoOnGeoType"));
501 const MEDCoupling::MEDFileMesh *m(0);
503 m=ms->getMeshWithName(meshName);
504 const MEDCoupling::MEDFileFastCellSupportComparator *fsp(_fsp);
505 if(!fsp || fsp->getNumberOfTS()==0)
507 std::vector< INTERP_KERNEL::NormalizedCellType > gts(fsp->getGeoTypesAt(0,m));
508 for(std::vector< INTERP_KERNEL::NormalizedCellType >::const_iterator it2=gts.begin();it2!=gts.end();it2++)
510 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(*it2));
511 std::string cmStr(cm.getRepr()); cmStr=cmStr.substr(5);//skip "NORM_"
512 sil->AddChild(root3,edge);
513 names.push_back(cmStr);
517 bool MEDFileFieldRepresentationLeaves::containId(int id) const
519 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
520 if((*it).getId()==id)
525 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
527 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
528 if((*it).getZeName()==name)
536 void MEDFileFieldRepresentationLeaves::dumpState(std::map<std::string,bool>& status) const
538 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
539 status[(*it).getZeName()]=(*it).getStatus();
542 bool MEDFileFieldRepresentationLeaves::isActivated() const
544 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
545 if((*it).getStatus())
550 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
552 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
554 os << " - " << (*it0).getZeName() << " (";
555 if((*it0).getStatus())
559 os << ")" << std::endl;
563 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
565 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
566 (*it).setStatus(true);
569 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
571 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
572 if((*it).getId()==id)
574 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
577 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
580 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
581 std::vector<double> ret;
582 std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
583 return tk.getTimeStepsRegardingPolicy(dtits,ret);
586 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
589 return _arrays[0]->getTimeSteps(ts);
593 return std::vector< std::pair<int,int> >();
597 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
599 std::ostringstream oss;
600 oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
601 std::vector<double> ret1;
602 std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
603 std::size_t sz(ret1.size());
604 for(std::size_t i=0;i<sz;i++)
606 oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
609 std::string tmp(oss.str());
610 if(tmp.size()>200 && i!=sz-1)
620 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshes *meshes, vtkDataSet *ds) const
623 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
624 MCAuto<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
625 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
626 if((*it).getStatus())
628 (*it).appendFields(tr,globs,mml,mst,ds);
629 (*it).appendELGAIfAny(ds);
633 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
635 const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
636 DataArrayDouble *coordsMC(0);
637 DataArrayByte *typesMC(0);
638 DataArrayInt *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
639 bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
640 MCAuto<DataArrayDouble> coordsSafe(coordsMC);
641 MCAuto<DataArrayByte> typesSafe(typesMC);
642 MCAuto<DataArrayInt> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
644 int nbOfCells(typesSafe->getNbOfElems());
645 vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
646 vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
647 cellTypes->SetArray(reinterpret_cast<unsigned char *>(typesSafe->getPointer()),nbOfCells,0,VTK_DATA_ARRAY_FREE); typesSafe->accessToMemArray().setSpecificDeallocator(0);
648 vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
649 cellLocations->SetVoidArray(cellLocationsSafe->getPointer(),nbOfCells,0,VTK_DATA_ARRAY_FREE); cellLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
650 vtkCellArray *cells(vtkCellArray::New());
651 vtkIdTypeArray *cells2(vtkIdTypeArray::New());
652 cells2->SetVoidArray(cellsSafe->getPointer(),cellsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); cellsSafe->accessToMemArray().setSpecificDeallocator(0);
653 cells->SetCells(nbOfCells,cells2);
655 if(faceLocationsMC!=0 && facesMC!=0)
657 vtkIdTypeArray *faces(vtkIdTypeArray::New());
658 faces->SetVoidArray(facesSafe->getPointer(),facesSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); facesSafe->accessToMemArray().setSpecificDeallocator(0);
659 vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
660 faceLocations->SetVoidArray(faceLocationsSafe->getPointer(),faceLocationsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); faceLocationsSafe->accessToMemArray().setSpecificDeallocator(0);
661 ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
662 faceLocations->Delete();
666 ret->SetCells(cellTypes,cellLocations,cells);
668 cellLocations->Delete();
670 vtkPoints *pts(vtkPoints::New());
671 vtkDoubleArray *pts2(vtkDoubleArray::New());
672 pts2->SetNumberOfComponents(3);
675 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);
676 coordsSafe->accessToMemArray().setSpecificDeallocator(0);
679 pts2->SetArray(coordsSafe->getPointer(),coordsSafe->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
688 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(MEDCoupling::MEDCMeshMultiLev *mm) const
690 const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
692 std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
693 vtkDoubleArray *vtkTmp(0);
694 vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
695 std::size_t dim(arrs.size());
697 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
698 int sizePerAxe[3]={1,1,1};
699 sizePerAxe[0]=arrs[0]->getNbOfElems();
701 sizePerAxe[1]=arrs[1]->getNbOfElems();
703 sizePerAxe[2]=arrs[2]->getNbOfElems();
704 ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
705 vtkTmp=vtkDoubleArray::New();
706 vtkTmp->SetNumberOfComponents(1);
708 vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
710 { vtkTmp->SetArray(arrs[0]->getPointer(),arrs[0]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[0]->accessToMemArray().setSpecificDeallocator(0); }
711 ret->SetXCoordinates(vtkTmp);
716 vtkTmp=vtkDoubleArray::New();
717 vtkTmp->SetNumberOfComponents(1);
719 vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
721 { vtkTmp->SetArray(arrs[1]->getPointer(),arrs[1]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[1]->accessToMemArray().setSpecificDeallocator(0); }
722 ret->SetYCoordinates(vtkTmp);
728 vtkTmp=vtkDoubleArray::New();
729 vtkTmp->SetNumberOfComponents(1);
731 vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
733 { vtkTmp->SetArray(arrs[2]->getPointer(),arrs[2]->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); arrs[2]->accessToMemArray().setSpecificDeallocator(0); }
734 ret->SetZCoordinates(vtkTmp);
741 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(MEDCoupling::MEDCurveLinearMeshMultiLev *mm) const
743 const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
744 int meshStr[3]={1,1,1};
745 DataArrayDouble *coords(0);
746 std::vector<int> nodeStrct;
748 mm->buildVTUArrays(coords,nodeStrct,isInternal);
749 std::size_t dim(nodeStrct.size());
751 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
752 meshStr[0]=nodeStrct[0];
754 meshStr[1]=nodeStrct[1];
756 meshStr[2]=nodeStrct[2];
757 vtkStructuredGrid *ret(vtkStructuredGrid::New());
758 ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
759 vtkDoubleArray *da(vtkDoubleArray::New());
760 da->SetNumberOfComponents(3);
761 if(coords->getNumberOfComponents()==3)
764 da->SetArray(coords->getPointer(),coords->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);//VTK has not the ownership of double * because MEDLoader main struct has it !
766 { da->SetArray(coords->getPointer(),coords->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); coords->accessToMemArray().setSpecificDeallocator(0); }
770 MCAuto<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
771 da->SetArray(coords2->getPointer(),coords2->getNbOfElems(),0,VTK_DATA_ARRAY_FREE);//let VTK deal with double *
772 coords2->accessToMemArray().setSpecificDeallocator(0);
775 vtkPoints *points=vtkPoints::New();
776 ret->SetPoints(points);
783 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes) const
785 const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
787 //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
788 MCAuto<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change
789 MCAuto<MEDMeshMultiLev> mml2(mml->prepare());
790 MEDMeshMultiLev *ptMML2(mml2);
793 MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
794 MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
795 MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
799 ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
803 ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
807 ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
810 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
811 _cached_ds=ret->NewInstance();
812 _cached_ds->ShallowCopy(ret);
816 ret=_cached_ds->NewInstance();
817 ret->ShallowCopy(_cached_ds);
820 appendFields(tr,globs,mml,meshes,ret);
821 // The arrays links to mesh
822 DataArrayInt *famCells(0),*numCells(0);
823 bool noCpyFamCells(false),noCpyNumCells(false);
824 ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
827 vtkIntArray *vtkTab(vtkIntArray::New());
828 vtkTab->SetNumberOfComponents(1);
829 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
831 vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
833 { vtkTab->SetArray(famCells->getPointer(),famCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famCells->accessToMemArray().setSpecificDeallocator(0); }
834 ret->GetCellData()->AddArray(vtkTab);
838 ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
841 vtkIntArray *vtkTab(vtkIntArray::New());
842 vtkTab->SetNumberOfComponents(1);
843 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
845 vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
847 { vtkTab->SetArray(numCells->getPointer(),numCells->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numCells->accessToMemArray().setSpecificDeallocator(0); }
848 ret->GetCellData()->AddArray(vtkTab);
852 // The arrays links to mesh
853 DataArrayInt *famNodes(0),*numNodes(0);
854 bool noCpyFamNodes(false),noCpyNumNodes(false);
855 ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
858 vtkIntArray *vtkTab(vtkIntArray::New());
859 vtkTab->SetNumberOfComponents(1);
860 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
862 vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
864 { vtkTab->SetArray(famNodes->getPointer(),famNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); famNodes->accessToMemArray().setSpecificDeallocator(0); }
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);
876 vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),1,VTK_DATA_ARRAY_FREE);
878 { vtkTab->SetArray(numNodes->getPointer(),numNodes->getNbOfElems(),0,VTK_DATA_ARRAY_FREE); numNodes->accessToMemArray().setSpecificDeallocator(0); }
879 ret->GetPointData()->AddArray(vtkTab);
886 //////////////////////
888 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
892 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
895 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
896 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
897 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
898 ret+=(*it2).getNumberOfArrays();
902 void MEDFileFieldRepresentationTree::assignIds() const
905 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
906 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
907 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
911 void MEDFileFieldRepresentationTree::computeFullNameInLeaves() const
913 std::size_t it0Cnt(0);
914 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
916 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
917 std::string tsName(oss.str());
918 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
920 std::string meshName((*it1)[0].getMeshName());
921 std::size_t it2Cnt(0);
922 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
924 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
925 std::string comSupStr(oss2.str());
926 (*it2).computeFullNameInLeaves(tsName,meshName,comSupStr);
932 void MEDFileFieldRepresentationTree::activateTheFirst() const
934 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
935 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
936 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
938 (*it2).activateAllArrays();
943 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
945 std::size_t it0Cnt(0);
946 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
948 vtkIdType InfoOnTSId(sil->AddChild(root,edge));
949 names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
951 vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
952 std::vector<double> ts;
953 std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
954 std::size_t nbOfTS(dtits.size());
955 std::ostringstream oss3; oss3 << nbOfTS;
956 names.push_back(oss3.str());
957 for(std::size_t i=0;i<nbOfTS;i++)
959 std::ostringstream oss4; oss4 << dtits[i].first;
960 vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
961 names.push_back(oss4.str());
962 std::ostringstream oss5; oss5 << dtits[i].second;
963 vtkIdType ItId(sil->AddChild(DtId,edge));
964 names.push_back(oss5.str());
965 std::ostringstream oss6; oss6 << ts[i];
966 sil->AddChild(ItId,edge);
967 names.push_back(oss6.str());
970 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
971 std::string tsName(oss.str());
972 vtkIdType typeId0(sil->AddChild(root,edge));
973 names.push_back(tsName);
974 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
976 std::string meshName((*it1)[0].getMeshName());
977 vtkIdType typeId1(sil->AddChild(typeId0,edge));
978 names.push_back(meshName);
979 std::size_t it2Cnt(0);
980 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
982 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
983 std::string comSupStr(oss2.str());
984 vtkIdType typeId2(sil->AddChild(typeId1,edge));
985 names.push_back(comSupStr);
986 (*it2).feedSIL(_ms,meshName,sil,typeId2,edge,names);
992 std::string MEDFileFieldRepresentationTree::getActiveMeshName() const
994 int dummy0(0),dummy1(0),dummy2(0);
995 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
996 return leaf.getMeshName();
999 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
1001 int dummy0(0),dummy1(0),dummy2(0);
1002 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
1003 std::string ret(leaf.getMeshName());
1006 for(;i<_ms->getNumberOfMeshes();i++)
1008 m=_ms->getMeshAtPos(i);
1009 if(m->getName()==ret)
1012 if(i==_ms->getNumberOfMeshes())
1013 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
1014 vtkIdType typeId0(sil->AddChild(root,edge));
1015 names.push_back(m->getName());
1017 vtkIdType typeId1(sil->AddChild(typeId0,edge));
1018 names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
1019 std::vector<std::string> grps(m->getGroupsNames());
1020 for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
1022 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1023 names.push_back(*it0);
1024 std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
1025 for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
1027 sil->AddChild(typeId2,edge);
1028 names.push_back((*it1).c_str());
1032 vtkIdType typeId11(sil->AddChild(typeId0,edge));
1033 names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
1034 std::vector<std::string> fams(m->getFamiliesNames());
1035 for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
1037 sil->AddChild(typeId11,edge);
1038 int famId(m->getFamilyId((*it00).c_str()));
1039 std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
1040 names.push_back(oss.str());
1045 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
1047 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1048 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1049 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1050 if((*it2).containId(id))
1051 return (*it2).getLeafArr(id);
1052 throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
1055 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
1057 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1058 return elt.getZeName();
1061 const char *MEDFileFieldRepresentationTree::getNameOfC(int id) const
1063 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1064 return elt.getZeNameC();
1067 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1069 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1070 return elt.getStatus();
1073 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1076 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1077 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1078 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1079 if((*it2).containZeName(name,ret))
1081 std::ostringstream msg; msg << "MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name \"" << name << "\" !";
1082 throw INTERP_KERNEL::Exception(msg.str().c_str());
1085 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1087 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1088 bool ret(elt.setStatus(status));//to be implemented
1092 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1095 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1096 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1097 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1098 ret=std::max(ret,(*it2).getNumberOfTS());
1106 void MEDFileFieldRepresentationTree::loadInMemory(MEDFileFields *fields, MEDFileMeshes *meshes)
1108 if(!((MEDCoupling::MEDFileFields *)fields))
1110 fields=BuildFieldFromMeshes(meshes);
1114 AppendFieldFromMeshes(meshes,fields);
1116 meshes->cartesianizeMe();
1117 fields->removeFieldsWithoutAnyTimeStep();
1118 std::vector<std::string> meshNames(meshes->getMeshesNames());
1119 std::vector< MCAuto<MEDFileFields> > fields_per_mesh(meshNames.size());
1120 for(std::size_t i=0;i<meshNames.size();i++)
1122 fields_per_mesh[i]=fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1124 std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1126 for(std::vector< MCAuto<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1128 for(int j=0;j<(*fields)->getNumberOfFields();j++)
1130 MCAuto<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1131 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1133 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::const_iterator it=tmp.begin();it!=tmp.end();it++)
1135 if(!(*it)->presenceOfMultiDiscPerGeoType())
1136 allFMTSLeavesToDisplaySafe.push_back(*it);
1138 {// The case of some parts of field have more than one discretization per geo type.
1139 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > subTmp((*it)->splitMultiDiscrPerGeoTypes());
1140 std::size_t it0Cnt(0);
1141 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::iterator it0=subTmp.begin();it0!=subTmp.end();it0++,it0Cnt++)//not const because setName
1143 std::ostringstream oss; oss << (*it0)->getName() << "_" << std::setfill('M') << std::setw(3) << it0Cnt;
1144 (*it0)->setName(oss.str());
1145 allFMTSLeavesToDisplaySafe.push_back(*it0);
1152 std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1153 for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1155 allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1157 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1158 // memory safety part
1159 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1160 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1162 allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1163 for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1165 allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1166 allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1169 // end of memory safety part
1170 // 1st : timesteps, 2nd : meshName, 3rd : common support
1171 this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1172 for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1174 std::vector< std::string > meshNamesLoc;
1175 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1176 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1178 std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1179 std::vector< std::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1180 if(it==meshNamesLoc.end())
1182 meshNamesLoc.push_back(meshName);
1183 splitByMeshName.resize(splitByMeshName.size()+1);
1184 splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1187 splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1189 _data_structure[i].resize(meshNamesLoc.size());
1190 for(std::size_t j=0;j<splitByMeshName.size();j++)
1192 std::vector< MCAuto<MEDFileFastCellSupportComparator> > fsp;
1193 std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1194 for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1195 sbmn[k]=splitByMeshName[j][k];
1196 //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1197 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,meshes->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1198 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1199 this->_data_structure[i][j].resize(commonSupSplit.size());
1200 for(std::size_t k=0;k<commonSupSplit.size();k++)
1202 commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1203 for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1205 commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1206 commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1209 for(std::size_t k=0;k<commonSupSplit.size();k++)
1210 this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1213 this->removeEmptyLeaves();
1215 this->computeFullNameInLeaves();
1219 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts)
1223 if((iPart==-1 && nbOfParts==-1) || (iPart==0 && nbOfParts==1))
1225 _ms=MEDFileMeshes::New(fileName);
1226 _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
1230 #ifdef MEDREADER_USE_MPI
1231 _ms=ParaMEDFileMeshes::New(iPart,nbOfParts,fileName);
1232 int nbMeshes(_ms->getNumberOfMeshes());
1233 for(int i=0;i<nbMeshes;i++)
1235 MEDCoupling::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1236 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1238 MCAuto<DataArrayInt> tmp3(tmp2->zipCoords());
1240 _fields=MEDFileFields::LoadPartOf(fileName,false,_ms);//false is important to not read the values
1242 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::loadMainStructureOfFile : request for iPart/nbOfParts=" << iPart << "/" << nbOfParts << " whereas Plugin not compiled with MPI !";
1243 throw INTERP_KERNEL::Exception(oss.str().c_str());
1249 MCAuto<MEDCoupling::SauvReader> sr(MEDCoupling::SauvReader::New(fileName));
1250 MCAuto<MEDCoupling::MEDFileData> mfd(sr->loadInMEDFileDS());
1251 _ms=mfd->getMeshes(); _ms->incrRef();
1252 int nbMeshes(_ms->getNumberOfMeshes());
1253 for(int i=0;i<nbMeshes;i++)
1255 MEDCoupling::MEDFileMesh *tmp(_ms->getMeshAtPos(i));
1256 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1258 tmp2->forceComputationOfParts();
1260 _fields=mfd->getFields();
1261 if((MEDCoupling::MEDFileFields *)_fields)
1265 loadInMemory(_fields,_ms);
1268 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1270 std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1271 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1273 std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1274 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1276 std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1277 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1279 newSD1.push_back(*it2);
1281 newSD0.push_back(newSD1);
1284 newSD.push_back(newSD0);
1288 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1290 if(compInfos.size()!=1)
1292 return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1295 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1297 return _data_structure[0][0][0].getMeshName();
1300 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1303 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1304 return leaf.getTimeSteps(tk);
1307 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const
1310 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1311 meshName=leaf.getMeshName();
1312 std::vector<double> ts(leaf.getTimeSteps(tk));
1313 std::size_t zeTimeId(0);
1316 std::vector<double> ts2(ts.size());
1317 std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1318 std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1319 zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1322 if(zeTimeId==(int)ts.size())
1323 zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1324 if(zeTimeId==(int)ts.size())
1325 {//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.
1326 //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1328 double valAttachedToPos(-std::numeric_limits<double>::max());
1329 for(std::size_t i=0;i<ts.size();i++)
1333 if(ts[i]>valAttachedToPos)
1336 valAttachedToPos=ts[i];
1341 {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1342 valAttachedToPos=std::numeric_limits<double>::max();
1343 for(std::size_t i=0;i<ts.size();i++)
1345 if(ts[i]<valAttachedToPos)
1348 valAttachedToPos=ts[i];
1353 std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1354 std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1355 oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1356 std::cerr << oss.str() << std::endl;
1360 tr=new MEDStdTimeReq((int)zeTimeId);
1362 tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
1363 vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms));
1368 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1370 int nbOfActivated(0);
1371 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1372 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1373 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1374 if((*it2).isActivated())
1376 if(nbOfActivated!=1)
1378 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1379 throw INTERP_KERNEL::Exception(oss.str().c_str());
1381 int i0(0),i1(0),i2(0);
1382 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1383 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1384 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1385 if((*it2).isActivated())
1387 lev0=i0; lev1=i1; lev2=i2;
1390 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1393 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1396 os << "#############################################" << std::endl;
1397 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1400 os << "TS" << i << std::endl;
1401 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1404 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1407 os << " " << (*it2).getMeshName() << std::endl;
1408 os << " Comp" << k << std::endl;
1409 (*it2).printMySelf(os);
1413 os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1416 std::map<std::string,bool> MEDFileFieldRepresentationTree::dumpState() const
1418 std::map<std::string,bool> ret;
1419 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1420 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1421 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1422 (*it2).dumpState(ret);
1426 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms, MEDCoupling::MEDFileFields *ret)
1429 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::AppendFieldFromMeshes : internal error ! NULL ret !");
1430 for(int i=0;i<ms->getNumberOfMeshes();i++)
1432 MEDFileMesh *mm(ms->getMeshAtPos(i));
1433 std::vector<int> levs(mm->getNonEmptyLevels());
1434 MEDCoupling::MCAuto<MEDCoupling::MEDFileField1TS> f1tsMultiLev(MEDCoupling::MEDFileField1TS::New());
1435 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1438 for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1440 std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1441 for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1443 MEDCoupling::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1444 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1446 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1447 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1450 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1451 f1tsMultiLev->setFieldNoProfileSBT(f);
1456 std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1457 if(levsExt.size()==levs.size()+1)
1459 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(1));
1460 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_NODES));
1462 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1463 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1464 arr->iota(); f->setArray(arr);
1465 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1466 f1tsMultiLev->setFieldNoProfileSBT(f);
1474 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(0));
1475 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1477 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1478 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1481 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1482 f1tsMultiLev->setFieldNoProfileSBT(f);
1485 MEDCoupling::MCAuto<MEDCoupling::MEDFileFieldMultiTS> fmtsMultiLev(MEDCoupling::MEDFileFieldMultiTS::New());
1486 fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1487 ret->pushField(fmtsMultiLev);
1491 std::string MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh(const std::string& meshName, const MEDCoupling::MEDFileFields *ret)
1493 const char KEY_STR_TO_AVOID_COLLIDE[]="MESH@";
1495 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh : internal error ! NULL ret !");
1496 std::vector<std::string> fieldNamesAlreadyExisting(ret->getFieldsNames());
1497 if(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),meshName)==fieldNamesAlreadyExisting.end())
1499 std::string tmpName(KEY_STR_TO_AVOID_COLLIDE); tmpName+=meshName;
1500 while(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),tmpName)!=fieldNamesAlreadyExisting.end())
1501 tmpName=std::string(KEY_STR_TO_AVOID_COLLIDE)+tmpName;
1505 MEDCoupling::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms)
1507 MEDCoupling::MCAuto<MEDCoupling::MEDFileFields> ret(MEDCoupling::MEDFileFields::New());
1508 AppendFieldFromMeshes(ms,ret);
1512 std::vector<std::string> MEDFileFieldRepresentationTree::SplitFieldNameIntoParts(const std::string& fullFieldName, char sep)
1514 std::vector<std::string> ret;
1516 while(pos!=std::string::npos)
1518 std::size_t curPos(fullFieldName.find_first_of(sep,pos));
1519 std::string elt(fullFieldName.substr(pos,curPos!=std::string::npos?curPos-pos:std::string::npos));
1521 pos=fullFieldName.find_first_not_of(sep,curPos);
1527 * Here the non regression tests.
1528 * const char inp0[]="";
1529 * const char exp0[]="";
1530 * const char inp1[]="field";
1531 * const char exp1[]="field";
1532 * const char inp2[]="_________";
1533 * const char exp2[]="_________";
1534 * const char inp3[]="field_p";
1535 * const char exp3[]="field_p";
1536 * const char inp4[]="field__p";
1537 * const char exp4[]="field_p";
1538 * const char inp5[]="field_p__";
1539 * const char exp5[]="field_p";
1540 * const char inp6[]="field_p_";
1541 * const char exp6[]="field_p";
1542 * const char inp7[]="field_____EDFGEG//sdkjf_____PP_______________";
1543 * const char exp7[]="field_EDFGEG//sdkjf_PP";
1544 * const char inp8[]="field_____EDFGEG//sdkjf_____PP";
1545 * const char exp8[]="field_EDFGEG//sdkjf_PP";
1546 * const char inp9[]="_field_____EDFGEG//sdkjf_____PP_______________";
1547 * const char exp9[]="field_EDFGEG//sdkjf_PP";
1548 * const char inp10[]="___field_____EDFGEG//sdkjf_____PP_______________";
1549 * const char exp10[]="field_EDFGEG//sdkjf_PP";
1551 std::string MEDFileFieldRepresentationTree::PostProcessFieldName(const std::string& fullFieldName)
1553 const char SEP('_');
1554 std::vector<std::string> v(SplitFieldNameIntoParts(fullFieldName,SEP));
1556 return fullFieldName;//should never happen
1560 return fullFieldName;
1564 std::string ret(v[0]);
1565 for(std::size_t i=1;i<v.size();i++)
1570 { ret+=SEP; ret+=v[i]; }
1576 return fullFieldName;
1582 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1586 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1591 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1593 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1595 throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1601 * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1602 * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1603 * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1605 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1607 std::size_t sz(ts.size());
1608 bool isInHumanRange(true);
1610 for(std::size_t i=0;i<sz;i++)
1613 if(ts[i]<=-1e299 || ts[i]>=1e299)
1614 isInHumanRange=false;
1617 return processedUsingPairOfIds(tsPairs);
1619 return processedUsingPairOfIds(tsPairs);
1620 _postprocessed_time=ts;
1621 return getPostProcessedTime();
1626 * idem than 0, except that ts is preaccumulated before invoking policy 0.
1628 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1630 std::size_t sz(ts.size());
1631 std::vector<double> ts2(sz);
1633 for(std::size_t i=0;i<sz;i++)
1638 return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1641 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1643 std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1647 void TimeKeeper::printSelf(std::ostream& oss) const
1649 std::size_t sz(_activated_ts.size());
1650 for(std::size_t i=0;i<sz;i++)
1652 oss << "(" << i << "," << _activated_ts[i].first << "), ";
1656 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1658 std::size_t sz(_activated_ts.size());
1659 std::vector<bool> ret(sz);
1660 for(std::size_t i=0;i<sz;i++)
1662 ret[i]=_activated_ts[i].first;
1667 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1669 std::size_t sz(tsPairs.size());
1670 std::set<int> s0,s1;
1671 for(std::size_t i=0;i<sz;i++)
1672 { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1675 _postprocessed_time.resize(sz);
1676 for(std::size_t i=0;i<sz;i++)
1677 _postprocessed_time[i]=(double)tsPairs[i].first;
1678 return getPostProcessedTime();
1682 _postprocessed_time.resize(sz);
1683 for(std::size_t i=0;i<sz;i++)
1684 _postprocessed_time[i]=(double)tsPairs[i].second;
1685 return getPostProcessedTime();
1687 //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1688 _postprocessed_time.resize(sz);
1689 for(std::size_t i=0;i<sz;i++)
1690 _postprocessed_time[i]=(double)i;
1691 return getPostProcessedTime();
1694 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1696 _activated_ts.resize(maxNumberOfTS);
1697 for(int i=0;i<maxNumberOfTS;i++)
1699 std::ostringstream oss; oss << "000" << i;
1700 _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());