1 // Copyright (C) 2010-2020 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"
30 #include "MEDCouplingMemArray.txx"
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 "vtkLongArray.h"
52 #include "vtkLongLongArray.h"
54 #include "vtkFloatArray.h"
55 #include "vtkCellArray.h"
56 #include "vtkPointData.h"
57 #include "vtkFieldData.h"
58 #include "vtkCellData.h"
60 #include "vtkMutableDirectedGraph.h"
62 using namespace MEDCoupling;
64 const char MEDFileFieldRepresentationLeavesArrays::ZE_SEP[]="@@][@@";
66 const char MEDFileFieldRepresentationLeavesArrays::TS_STR[]="TS";
68 const char MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR[]="ComSup";
70 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME[]="FamilyIdCell";
72 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME[]="NumIdCell";
74 const char MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME[]="FamilyIdNode";
76 const char MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME[]="NumIdNode";
78 const char MEDFileFieldRepresentationLeavesArrays::GLOBAL_NODE_ID_NAME[]="GlobalNodeIds";// WARNING DO NOT CHANGE IT BEFORE HAVING CHECKED IN PV SOURCES !
80 const char MEDFileFieldRepresentationTree::ROOT_OF_GRPS_IN_TREE[]="zeGrps";
82 const char MEDFileFieldRepresentationTree::ROOT_OF_FAM_IDS_IN_TREE[]="zeFamIds";
84 const char MEDFileFieldRepresentationTree::COMPO_STR_TO_LOCATE_MESH_DA[]="-@?|*_";
87 vtkIdTypeArray *ELGACmp::findOrCreate(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDataArray *vtkd, vtkDataSet *ds, bool& isNew, ExportedTinyInfo *internalInfo) const
89 vtkIdTypeArray *try0(isExisting(locsReallyUsed,vtkd));
98 return createNew<T>(globs,locsReallyUsed,vtkd,ds,internalInfo);
102 vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUsed, vtkDataArray *vtkd) const
104 std::vector< std::vector<std::string> >::iterator it(std::find(_loc_names.begin(),_loc_names.end(),locsReallyUsed));
105 if(it==_loc_names.end())
107 std::size_t pos(std::distance(_loc_names.begin(),it));
108 vtkIdTypeArray *ret(_elgas[pos]);
109 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
110 for(std::vector<std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it=_defs[pos].begin();it!=_defs[pos].end();it++)
112 key->Set(vtkd->GetInformation(),(*it).first,(*it).second);
114 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),ret->GetName());
119 vtkIdTypeArray *ELGACmp::createNew(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDataArray *vtkd, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
121 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
122 std::vector< std::vector<std::string> > locNames(_loc_names);
123 std::vector<vtkIdTypeArray *> elgas(_elgas);
124 std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > defs;
126 std::vector< std::vector<std::string> >::const_iterator it(std::find(locNames.begin(),locNames.end(),locsReallyUsed));
127 if(it!=locNames.end())
128 throw INTERP_KERNEL::Exception("ELGACmp::createNew : Method is expected to be called after isExisting call ! Entry already exists !");
129 locNames.push_back(locsReallyUsed);
130 vtkIdTypeArray *elga(vtkIdTypeArray::New());
131 elga->SetNumberOfComponents(1);
132 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
133 std::map<unsigned char,int> m;
134 for(std::vector<std::string>::const_iterator it=locsReallyUsed.begin();it!=locsReallyUsed.end();it++)
136 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
137 const MEDFileFieldLoc& loc(globs->getLocalization((*it).c_str()));
138 INTERP_KERNEL::NormalizedCellType ct(loc.getGeoType());
139 unsigned char vtkType(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE[ct]);
140 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
141 int nbGaussPt(loc.getNbOfGaussPtPerCell()),nbPtsPerCell((int)cm.getNumberOfNodes()),dimLoc(loc.getDimension());
142 // 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.
143 std::vector<double> gsCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getGaussCoords()));
144 std::vector<double> refCoods2(INTERP_KERNEL::GaussInfo::NormalizeCoordinatesIfNecessary(ct,dimLoc,loc.getRefCoords()));
146 internalInfo->pushGaussAdditionnalInfo(vtkType,dimLoc,refCoods2,gsCoods2);
147 double *shape(new double[nbPtsPerCell*nbGaussPt]);
148 INTERP_KERNEL::GaussInfo calculator(ct,gsCoods2,nbGaussPt,refCoods2,nbPtsPerCell);
149 calculator.initLocalInfo();
150 const std::vector<double>& wgths(loc.getGaussWeights());
151 for(int i=0;i<nbGaussPt;i++)
153 const double *pt0(calculator.getFunctionValues(i));
154 if(ct!=INTERP_KERNEL::NORM_HEXA27)
155 std::copy(pt0,pt0+nbPtsPerCell,shape+nbPtsPerCell*i);
158 for(int j=0;j<27;j++)
159 shape[nbPtsPerCell*i+j]=pt0[MEDMeshMultiLev::HEXA27_PERM_ARRAY[j]];
162 m[vtkType]=nbGaussPt;
163 def->Initialize(vtkType,nbPtsPerCell,nbGaussPt,shape,const_cast<double *>(&wgths[0]));
165 key->Set(elga->GetInformation(),def,vtkType);
166 key->Set(vtkd->GetInformation(),def,vtkType);
167 defs.push_back(std::pair< vtkQuadratureSchemeDefinition *, unsigned char >(def,vtkType));
170 vtkIdType ncell(ds->GetNumberOfCells());
171 vtkIdType *pt(new vtkIdType[ncell]),offset(0);
172 for(vtkIdType cellId=0;cellId<ncell;cellId++)
174 vtkCell *cell(ds->GetCell(cellId));
175 vtkIdType delta(m[(unsigned char)cell->GetCellType()]);
179 elga->GetInformation()->Set(MEDUtilities::ELGA(),1);
180 elga->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
181 std::ostringstream oss; oss << "ELGA" << "@" << _loc_names.size();
182 std::string ossStr(oss.str());
183 elga->SetName(ossStr.c_str());
184 elga->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
185 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elga->GetName());
186 elgas.push_back(elga);
190 _defs.push_back(defs);
194 void ELGACmp::appendELGAIfAny(vtkDataSet *ds) const
196 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
197 ds->GetCellData()->AddArray(*it);
202 for(std::vector<vtkIdTypeArray *>::const_iterator it=_elgas.begin();it!=_elgas.end();it++)
204 for(std::vector< std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > >::const_iterator it0=_defs.begin();it0!=_defs.end();it0++)
205 for(std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
206 (*it1).first->Delete();
212 class MEDFileVTKTraits
215 typedef void VtkType;
220 class MEDFileVTKTraits<int>
223 typedef vtkIntArray VtkType;
224 typedef MEDCoupling::DataArrayInt MCType;
229 class MEDFileVTKTraits<long long>
231 class MEDFileVTKTraits<long>
236 typedef vtkLongLongArray VtkType;
238 typedef vtkLongArray VtkType;
240 typedef MEDCoupling::DataArrayInt64 MCType;
244 class MEDFileVTKTraits<float>
247 typedef vtkFloatArray VtkType;
248 typedef MEDCoupling::DataArrayFloat MCType;
252 class MEDFileVTKTraits<double>
255 typedef vtkDoubleArray VtkType;
256 typedef MEDCoupling::DataArrayDouble MCType;
259 typedef typename MEDFileVTKTraits<mcIdType>::VtkType vtkMCIdTypeArray;
263 void AssignDataPointerToVTK(typename MEDFileVTKTraits<T>::VtkType *vtkTab, typename MEDFileVTKTraits<T>::MCType *mcTab, bool noCpyNumNodes)
266 vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),1,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE);
268 { vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),0,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE); mcTab->accessToMemArray().setSpecificDeallocator(0); }
271 // here copy is always assumed.
272 template<class VTKT, class MCT>
273 void AssignDataPointerOther(VTKT *vtkTab, MCT *mcTab, int nbElems)
275 typedef typename VTKT::ValueType VTKType;
276 if ( sizeof( VTKType ) == sizeof( typename MCT::Type ))
278 vtkTab->SetVoidArray(reinterpret_cast<unsigned char *>(mcTab->getPointer()),nbElems,0,VTKT::VTK_DATA_ARRAY_FREE);
279 mcTab->accessToMemArray().setSpecificDeallocator(0);
283 VTKType* newArray = new VTKType[ nbElems ];
284 std::copy( mcTab->begin(), mcTab->begin() + nbElems, newArray );
285 vtkTab->SetVoidArray(reinterpret_cast<unsigned char *>(newArray),nbElems,0,VTKT::VTK_DATA_ARRAY_DELETE);
290 void AssignToFieldData(DataArray *vPtr, const MEDTimeReq *tr, vtkFieldData *att, const std::string& crudeName, bool noCpyNumNodes,
291 const std::vector<TypeOfField>& discs, const ELGACmp& elgaCmp, const MEDCoupling::MEDFileFieldGlobsReal *globs,
292 MEDFileAnyTypeField1TS *f1ts, vtkDataSet *ds, ExportedTinyInfo *internalInfo)
294 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
295 typename MEDFileVTKTraits<T>::MCType *vi(static_cast<typename MEDFileVTKTraits<T>::MCType *>(vPtr));
296 typename MEDFileVTKTraits<T>::VtkType *vtkd(MEDFileVTKTraits<T>::VtkType::New());
297 vtkd->SetNumberOfComponents((int)vi->getNumberOfComponents());
298 for(unsigned int i=0;i<vi->getNumberOfComponents();i++)
299 vtkd->SetComponentName(i,vi->getVarOnComponent(i).c_str());
300 AssignDataPointerToVTK<T>(vtkd,vi,noCpyNumNodes);
301 std::string name(tr->buildName(crudeName));
302 vtkd->SetName(name.c_str());
305 if(discs[0]==ON_GAUSS_PT)
308 elgaCmp.findOrCreate<T>(globs,f1ts->getLocsReallyUsed(),vtkd,ds,tmp,internalInfo);
310 if(discs[0]==ON_GAUSS_NE)
312 vtkIdTypeArray *elno(vtkIdTypeArray::New());
313 elno->SetNumberOfComponents(1);
314 vtkIdType ncell(ds->GetNumberOfCells());
315 vtkIdType *pt(new vtkIdType[ncell]),offset(0);
316 std::set<int> cellTypes;
317 for(vtkIdType cellId=0;cellId<ncell;cellId++)
319 vtkCell *cell(ds->GetCell(cellId));
320 vtkIdType delta(cell->GetNumberOfPoints());
321 cellTypes.insert(cell->GetCellType());
325 elno->GetInformation()->Set(MEDUtilities::ELNO(),1);
326 elno->SetVoidArray(pt,ncell,0,VTK_DATA_ARRAY_DELETE);
327 std::string nameElno("ELNO"); nameElno+="@"; nameElno+=name;
328 elno->SetName(nameElno.c_str());
329 ds->GetCellData()->AddArray(elno);
330 vtkd->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),elno->GetName());
331 elno->GetInformation()->Set(vtkAbstractArray::GUI_HIDE(),1);
333 vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
334 for(std::set<int>::const_iterator it=cellTypes.begin();it!=cellTypes.end();it++)
336 const unsigned char *pos(std::find(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,*it));
337 if(pos==MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
339 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)std::distance(MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos));
340 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(ct));
341 int nbGaussPt(cm.getNumberOfNodes()),dim(cm.getDimension());
342 vtkQuadratureSchemeDefinition *def(vtkQuadratureSchemeDefinition::New());
343 double *shape(new double[nbGaussPt*nbGaussPt]);
345 const double *gsCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));//GetLocsFromGeometricType
346 const double *refCoords(MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(ct,dummy));
347 const double *weights(MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(ct,dummy));
348 std::vector<double> gsCoords2(gsCoords,gsCoords+nbGaussPt*dim),refCoords2(refCoords,refCoords+nbGaussPt*dim);
349 INTERP_KERNEL::GaussInfo calculator(ct,gsCoords2,nbGaussPt,refCoords2,nbGaussPt);
350 calculator.initLocalInfo();
351 for(int i=0;i<nbGaussPt;i++)
353 const double *pt0(calculator.getFunctionValues(i));
354 std::copy(pt0,pt0+nbGaussPt,shape+nbGaussPt*i);
356 def->Initialize(*it,nbGaussPt,nbGaussPt,shape,const_cast<double *>(weights));
358 key->Set(elno->GetInformation(),def,*it);
359 key->Set(vtkd->GetInformation(),def,*it);
369 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays():_id(-1)
373 MEDFileFieldRepresentationLeavesArrays::MEDFileFieldRepresentationLeavesArrays(const MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>& arr):MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>(arr),_activated(false),_id(-1)
375 std::vector< std::vector<MEDCoupling::TypeOfField> > typs((operator->())->getTypesOfFieldAvailable());
377 throw INTERP_KERNEL::Exception("There is a big internal problem in MEDLoader ! The field time spitting has failed ! A CRASH will occur soon !");
378 if(typs[0].size()!=1)
379 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 !");
380 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDiscretization> fd(MEDCouplingFieldDiscretization::New(typs[0][0]));
381 std::ostringstream oss2; oss2 << (operator->())->getName() << ZE_SEP << fd->getRepr();
385 MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeavesArrays::operator=(const MEDFileFieldRepresentationLeavesArrays& other)
387 MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS>::operator=(other);
390 _ze_name=other._ze_name;
391 _ze_full_name.clear();
395 void MEDFileFieldRepresentationLeavesArrays::setId(int& id) const
400 int MEDFileFieldRepresentationLeavesArrays::getId() const
405 std::string MEDFileFieldRepresentationLeavesArrays::getZeName() const
407 return _ze_full_name;
410 const char *MEDFileFieldRepresentationLeavesArrays::getZeNameC() const
412 return _ze_full_name.c_str();
415 void MEDFileFieldRepresentationLeavesArrays::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
417 vtkIdType refId(sil->AddChild(root,edge));
418 names.push_back(_ze_name);
420 if(MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(((operator->())->getInfo())))
422 sil->AddChild(refId,edge);
423 names.push_back(std::string());
427 void MEDFileFieldRepresentationLeavesArrays::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
429 std::ostringstream oss3; oss3 << tsName << "/" << meshName << "/" << comSupStr << "/" << _ze_name;
430 _ze_full_name=oss3.str();
433 bool MEDFileFieldRepresentationLeavesArrays::getStatus() const
438 bool MEDFileFieldRepresentationLeavesArrays::setStatus(bool status) const
440 bool ret(_activated!=status);
445 void MEDFileFieldRepresentationLeavesArrays::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshStruct *mst, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
447 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
448 tr->setNumberOfTS((operator->())->getNumberOfTS());
450 for(int timeStepId=0;timeStepId<tr->size();timeStepId++,++(*tr))
452 MCAuto<MEDFileAnyTypeField1TS> f1ts((operator->())->getTimeStepAtPos(tr->getCurrent()));
453 MEDFileAnyTypeField1TS *f1tsPtr(f1ts);
454 MEDFileField1TS *f1tsPtrDbl(dynamic_cast<MEDFileField1TS *>(f1tsPtr));
455 MEDFileInt32Field1TS *f1tsPtrInt(dynamic_cast<MEDFileInt32Field1TS *>(f1tsPtr));
456 MEDFileInt64Field1TS *f1tsPtrInt64(dynamic_cast<MEDFileInt64Field1TS *>(f1tsPtr));
457 MEDFileFloatField1TS *f1tsPtrFloat(dynamic_cast<MEDFileFloatField1TS *>(f1tsPtr));
458 DataArray *crudeArr(0),*postProcessedArr(0);
460 crudeArr=f1tsPtrDbl->getUndergroundDataArray();
462 crudeArr=f1tsPtrInt->getUndergroundDataArray();
463 else if(f1tsPtrInt64)
464 crudeArr=f1tsPtrInt64->getUndergroundDataArray();
465 else if(f1tsPtrFloat)
466 crudeArr=f1tsPtrFloat->getUndergroundDataArray();
468 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64, FLOAT32 and INT32 fields are dealt for the moment !");
469 MEDFileField1TSStructItem fsst(MEDFileField1TSStructItem::BuildItemFrom(f1ts,mst));
470 f1ts->loadArraysIfNecessary();
471 MCAuto<DataArray> v(mml->buildDataArray(fsst,globs,crudeArr));
474 std::vector<TypeOfField> discs(f1ts->getTypesOfFieldAvailable());
476 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : internal error ! Number of spatial discretizations must be equal to one !");
477 vtkFieldData *att(0);
482 att=ds->GetCellData();
487 att=ds->GetPointData();
492 att=ds->GetFieldData();
497 att=ds->GetFieldData();
501 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only CELL and NODE, GAUSS_NE and GAUSS fields are available for the moment !");
505 AssignToFieldData<double>(v,tr,att,f1ts->getName(),postProcessedArr==crudeArr,discs,_elga_cmp,globs,f1ts,ds,internalInfo);
509 AssignToFieldData<int>(v,tr,att,f1ts->getName(),postProcessedArr==crudeArr,discs,_elga_cmp,globs,f1ts,ds,internalInfo);
511 else if(f1tsPtrFloat)
513 AssignToFieldData<float>(v,tr,att,f1ts->getName(),postProcessedArr==crudeArr,discs,_elga_cmp,globs,f1ts,ds,internalInfo);
515 else if(f1tsPtrInt64)
517 AssignToFieldData<Int64>(v,tr,att,f1ts->getName(),postProcessedArr==crudeArr,discs,_elga_cmp,globs,f1ts,ds,internalInfo);
520 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeavesArrays::appendFields : only FLOAT64 and INT32 fields are dealt for the moment ! Internal Error !");
524 void MEDFileFieldRepresentationLeavesArrays::appendELGAIfAny(vtkDataSet *ds) const
526 _elga_cmp.appendELGAIfAny(ds);
531 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves():_cached_ds(0)
535 MEDFileFieldRepresentationLeaves::MEDFileFieldRepresentationLeaves(const std::vector< MEDCoupling::MCAuto<MEDCoupling::MEDFileAnyTypeFieldMultiTS> >& arr,
536 const MEDCoupling::MCAuto<MEDCoupling::MEDFileFastCellSupportComparator>& fsp):_arrays(arr.size()),_fsp(fsp),_cached_ds(0)
538 for(std::size_t i=0;i<arr.size();i++)
539 _arrays[i]=MEDFileFieldRepresentationLeavesArrays(arr[i]);
542 MEDFileFieldRepresentationLeaves::~MEDFileFieldRepresentationLeaves()
545 _cached_ds->Delete();
548 bool MEDFileFieldRepresentationLeaves::empty() const
550 const MEDFileFastCellSupportComparator *fcscp(_fsp);
551 return fcscp==0 || _arrays.empty();
554 void MEDFileFieldRepresentationLeaves::setId(int& id) const
556 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
560 std::string MEDFileFieldRepresentationLeaves::getMeshName() const
562 return _arrays[0]->getMeshName();
565 int MEDFileFieldRepresentationLeaves::getNumberOfArrays() const
567 return (int)_arrays.size();
570 int MEDFileFieldRepresentationLeaves::getNumberOfTS() const
572 return _arrays[0]->getNumberOfTS();
575 void MEDFileFieldRepresentationLeaves::computeFullNameInLeaves(const std::string& tsName, const std::string& meshName, const std::string& comSupStr) const
577 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
578 (*it).computeFullNameInLeaves(tsName,meshName,comSupStr);
582 * \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.
584 void MEDFileFieldRepresentationLeaves::feedSIL(const MEDCoupling::MEDFileMeshes *ms, const std::string& meshName, vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
586 vtkIdType root2(sil->AddChild(root,edge));
587 names.push_back(std::string("Arrs"));
588 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
589 (*it).feedSIL(sil,root2,edge,names);
591 vtkIdType root3(sil->AddChild(root,edge));
592 names.push_back(std::string("InfoOnGeoType"));
593 const MEDCoupling::MEDFileMesh *m(0);
595 m=ms->getMeshWithName(meshName);
596 const MEDCoupling::MEDFileFastCellSupportComparator *fsp(_fsp);
597 if(!fsp || fsp->getNumberOfTS()==0)
599 std::vector< INTERP_KERNEL::NormalizedCellType > gts(fsp->getGeoTypesAt(0,m));
600 for(std::vector< INTERP_KERNEL::NormalizedCellType >::const_iterator it2=gts.begin();it2!=gts.end();it2++)
602 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(*it2));
603 std::string cmStr(cm.getRepr()); cmStr=cmStr.substr(5);//skip "NORM_"
604 sil->AddChild(root3,edge);
605 names.push_back(cmStr);
609 bool MEDFileFieldRepresentationLeaves::containId(int id) const
611 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
612 if((*it).getId()==id)
617 bool MEDFileFieldRepresentationLeaves::containZeName(const char *name, int& id) const
619 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
620 if((*it).getZeName()==name)
628 void MEDFileFieldRepresentationLeaves::dumpState(std::map<std::string,bool>& status) const
630 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
631 status[(*it).getZeName()]=(*it).getStatus();
634 bool MEDFileFieldRepresentationLeaves::isActivated() const
636 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
637 if((*it).getStatus())
642 void MEDFileFieldRepresentationLeaves::printMySelf(std::ostream& os) const
644 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it0=_arrays.begin();it0!=_arrays.end();it0++)
646 os << " - " << (*it0).getZeName() << " (";
647 if((*it0).getStatus())
651 os << ")" << std::endl;
655 void MEDFileFieldRepresentationLeaves::activateAllArrays() const
657 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
658 (*it).setStatus(true);
661 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationLeaves::getLeafArr(int id) const
663 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
664 if((*it).getId()==id)
666 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getLeafArr ! No such id !");
669 std::vector<double> MEDFileFieldRepresentationLeaves::getTimeSteps(const TimeKeeper& tk) const
672 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::getTimeSteps : the array size must be at least of size one !");
673 std::vector<double> ret;
674 std::vector< std::pair<int,int> > dtits(_arrays[0]->getTimeSteps(ret));
675 return tk.getTimeStepsRegardingPolicy(dtits,ret);
678 std::vector< std::pair<int,int> > MEDFileFieldRepresentationLeaves::getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const
681 return _arrays[0]->getTimeSteps(ts);
685 return std::vector< std::pair<int,int> >();
689 std::string MEDFileFieldRepresentationLeaves::getHumanReadableOverviewOfTS() const
691 std::ostringstream oss;
692 oss << _arrays[0]->getNumberOfTS() << " time steps [" << _arrays[0]->getDtUnit() << "]\n(";
693 std::vector<double> ret1;
694 std::vector< std::pair<int,int> > ret2(getTimeStepsInCoarseMEDFileFormat(ret1));
695 std::size_t sz(ret1.size());
696 for(std::size_t i=0;i<sz;i++)
698 oss << ret1[i] << " (" << ret2[i].first << "," << ret2[i].second << ")";
701 std::string tmp(oss.str());
702 if(tmp.size()>200 && i!=sz-1)
712 void MEDFileFieldRepresentationLeaves::appendFields(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDMeshMultiLev *mml, const MEDCoupling::MEDFileMeshes *meshes, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
715 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::appendFields : internal error !");
716 MCAuto<MEDFileMeshStruct> mst(MEDFileMeshStruct::New(meshes->getMeshWithName(_arrays[0]->getMeshName().c_str())));
717 for(std::vector<MEDFileFieldRepresentationLeavesArrays>::const_iterator it=_arrays.begin();it!=_arrays.end();it++)
718 if((*it).getStatus())
720 (*it).appendFields(tr,globs,mml,mst,ds,internalInfo);
721 (*it).appendELGAIfAny(ds);
725 vtkUnstructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationUnstructured(MEDUMeshMultiLev *mm) const
727 DataArrayDouble *coordsMC(0);
728 DataArrayByte *typesMC(0);
729 DataArrayIdType *cellLocationsMC(0),*cellsMC(0),*faceLocationsMC(0),*facesMC(0);
730 bool statusOfCoords(mm->buildVTUArrays(coordsMC,typesMC,cellLocationsMC,cellsMC,faceLocationsMC,facesMC));
731 MCAuto<DataArrayDouble> coordsSafe(coordsMC);
732 MCAuto<DataArrayByte> typesSafe(typesMC);
733 MCAuto<DataArrayIdType> cellLocationsSafe(cellLocationsMC),cellsSafe(cellsMC),faceLocationsSafe(faceLocationsMC),facesSafe(facesMC);
735 int nbOfCells(typesSafe->getNbOfElems());
736 vtkUnstructuredGrid *ret(vtkUnstructuredGrid::New());
737 vtkUnsignedCharArray *cellTypes(vtkUnsignedCharArray::New());
738 AssignDataPointerOther<vtkUnsignedCharArray,DataArrayByte>(cellTypes,typesSafe,nbOfCells);
739 vtkIdTypeArray *cellLocations(vtkIdTypeArray::New());
740 AssignDataPointerOther<vtkIdTypeArray,DataArrayIdType>(cellLocations,cellLocationsSafe,nbOfCells);
741 vtkCellArray *cells(vtkCellArray::New());
742 vtkIdTypeArray *cells2(vtkIdTypeArray::New());
743 AssignDataPointerOther<vtkIdTypeArray,DataArrayIdType>(cells2,cellsSafe,cellsSafe->getNbOfElems());
744 cells->SetCells(nbOfCells,cells2);
746 if(faceLocationsMC!=0 && facesMC!=0)
748 vtkIdTypeArray *faces(vtkIdTypeArray::New());
749 AssignDataPointerOther<vtkIdTypeArray,DataArrayIdType>(faces,facesSafe,facesSafe->getNbOfElems());
750 vtkIdTypeArray *faceLocations(vtkIdTypeArray::New());
751 AssignDataPointerOther<vtkIdTypeArray,DataArrayIdType>(faceLocations,faceLocationsSafe,faceLocationsSafe->getNbOfElems());
752 ret->SetCells(cellTypes,cellLocations,cells,faceLocations,faces);
753 faceLocations->Delete();
757 ret->SetCells(cellTypes,cellLocations,cells);
759 cellLocations->Delete();
761 vtkPoints *pts(vtkPoints::New());
762 vtkDoubleArray *pts2(vtkDoubleArray::New());
763 pts2->SetNumberOfComponents(3);
764 AssignDataPointerToVTK<double>(pts2,coordsSafe,statusOfCoords);
773 vtkRectilinearGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCartesian(MEDCoupling::MEDCMeshMultiLev *mm) const
776 std::vector< DataArrayDouble * > arrs(mm->buildVTUArrays(isInternal));
777 vtkDoubleArray *vtkTmp(0);
778 vtkRectilinearGrid *ret(vtkRectilinearGrid::New());
779 std::size_t dim(arrs.size());
781 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCartesian : dimension must be in [1,3] !");
782 int sizePerAxe[3]={1,1,1};
783 sizePerAxe[0]=arrs[0]->getNbOfElems();
785 sizePerAxe[1]=arrs[1]->getNbOfElems();
787 sizePerAxe[2]=arrs[2]->getNbOfElems();
788 ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
789 vtkTmp=vtkDoubleArray::New();
790 vtkTmp->SetNumberOfComponents(1);
791 AssignDataPointerToVTK<double>(vtkTmp,arrs[0],isInternal);
792 ret->SetXCoordinates(vtkTmp);
797 vtkTmp=vtkDoubleArray::New();
798 vtkTmp->SetNumberOfComponents(1);
799 AssignDataPointerToVTK<double>(vtkTmp,arrs[1],isInternal);
800 ret->SetYCoordinates(vtkTmp);
806 vtkTmp=vtkDoubleArray::New();
807 vtkTmp->SetNumberOfComponents(1);
808 AssignDataPointerToVTK<double>(vtkTmp,arrs[2],isInternal);
809 ret->SetZCoordinates(vtkTmp);
816 vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolationCurveLinear(MEDCoupling::MEDCurveLinearMeshMultiLev *mm) const
818 int meshStr[3]={1,1,1};
819 DataArrayDouble *coords(0);
820 std::vector<mcIdType> nodeStrct;
822 mm->buildVTUArrays(coords,nodeStrct,isInternal);
823 std::size_t dim(nodeStrct.size());
825 throw INTERP_KERNEL::Exception("buildVTKInstanceNoTimeInterpolationCurveLinear : dimension must be in [1,3] !");
826 meshStr[0]=nodeStrct[0];
828 meshStr[1]=nodeStrct[1];
830 meshStr[2]=nodeStrct[2];
831 vtkStructuredGrid *ret(vtkStructuredGrid::New());
832 ret->SetDimensions(meshStr[0],meshStr[1],meshStr[2]);
833 vtkDoubleArray *da(vtkDoubleArray::New());
834 da->SetNumberOfComponents(3);
835 if(coords->getNumberOfComponents()==3)
836 AssignDataPointerToVTK<double>(da,coords,isInternal);//if isIntenal==True VTK has not the ownership of double * because MEDLoader main struct has it !
839 MCAuto<DataArrayDouble> coords2(coords->changeNbOfComponents(3,0.));
840 AssignDataPointerToVTK<double>(da,coords2,false);//let VTK deal with double *
843 vtkPoints *points=vtkPoints::New();
844 ret->SetPoints(points);
851 vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, ExportedTinyInfo *internalInfo) const
854 //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
855 MCAuto<MEDMeshMultiLev> mml(_fsp->buildFromScratchDataSetSupport(0,globs));//0=timestep Id. Make the hypothesis that support does not change
856 MCAuto<MEDMeshMultiLev> mml2(mml->prepare());
857 MEDMeshMultiLev *ptMML2(mml2);
860 MEDUMeshMultiLev *ptUMML2(dynamic_cast<MEDUMeshMultiLev *>(ptMML2));
861 MEDCMeshMultiLev *ptCMML2(dynamic_cast<MEDCMeshMultiLev *>(ptMML2));
862 MEDCurveLinearMeshMultiLev *ptCLMML2(dynamic_cast<MEDCurveLinearMeshMultiLev *>(ptMML2));
866 ret=buildVTKInstanceNoTimeInterpolationUnstructured(ptUMML2);
870 ret=buildVTKInstanceNoTimeInterpolationCartesian(ptCMML2);
874 ret=buildVTKInstanceNoTimeInterpolationCurveLinear(ptCLMML2);
877 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation : unrecognized mesh ! Supported for the moment unstructured, cartesian, curvelinear !");
878 _cached_ds=ret->NewInstance();
879 _cached_ds->ShallowCopy(ret);
883 ret=_cached_ds->NewInstance();
884 ret->ShallowCopy(_cached_ds);
887 appendFields(tr,globs,mml,meshes,ret,internalInfo);
888 // The arrays links to mesh
889 DataArrayIdType *famCells(0),*numCells(0);
890 bool noCpyFamCells(false),noCpyNumCells(false);
891 ptMML2->retrieveFamilyIdsOnCells(famCells,noCpyFamCells);
894 vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
895 vtkTab->SetNumberOfComponents(1);
896 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME);
897 AssignDataPointerToVTK<mcIdType>(vtkTab,famCells,noCpyFamCells);
898 ret->GetCellData()->AddArray(vtkTab);
902 ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
905 vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
906 vtkTab->SetNumberOfComponents(1);
907 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_CELL_NAME);
908 AssignDataPointerToVTK<mcIdType>(vtkTab,numCells,noCpyNumCells);
909 ret->GetCellData()->AddArray(vtkTab);
913 // The arrays links to mesh
914 DataArrayIdType *famNodes(0),*numNodes(0);
915 bool noCpyFamNodes(false),noCpyNumNodes(false);
916 ptMML2->retrieveFamilyIdsOnNodes(famNodes,noCpyFamNodes);
919 vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
920 vtkTab->SetNumberOfComponents(1);
921 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
922 AssignDataPointerToVTK<mcIdType>(vtkTab,famNodes,noCpyFamNodes);
923 ret->GetPointData()->AddArray(vtkTab);
927 ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
930 vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
931 vtkTab->SetNumberOfComponents(1);
932 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::NUM_ID_NODE_NAME);
933 AssignDataPointerToVTK<mcIdType>(vtkTab,numNodes,noCpyNumNodes);
934 ret->GetPointData()->AddArray(vtkTab);
938 // Global Node Ids if any ! (In // mode)
939 DataArrayIdType *gni(ptMML2->retrieveGlobalNodeIdsIfAny());
942 vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
943 vtkTab->SetNumberOfComponents(1);
944 vtkTab->SetName(MEDFileFieldRepresentationLeavesArrays::GLOBAL_NODE_ID_NAME);
945 AssignDataPointerToVTK<mcIdType>(vtkTab,gni,false);
946 ret->GetPointData()->AddArray(vtkTab);
953 //////////////////////
955 MEDFileFieldRepresentationTree::MEDFileFieldRepresentationTree()
959 int MEDFileFieldRepresentationTree::getNumberOfLeavesArrays() const
962 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
963 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
964 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
965 ret+=(*it2).getNumberOfArrays();
969 void MEDFileFieldRepresentationTree::assignIds() const
972 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
973 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
974 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
978 void MEDFileFieldRepresentationTree::computeFullNameInLeaves() const
980 std::size_t it0Cnt(0);
981 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
983 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
984 std::string tsName(oss.str());
985 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
987 std::string meshName((*it1)[0].getMeshName());
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 (*it2).computeFullNameInLeaves(tsName,meshName,comSupStr);
999 void MEDFileFieldRepresentationTree::activateTheFirst() const
1001 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1002 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1003 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1005 (*it2).activateAllArrays();
1010 void MEDFileFieldRepresentationTree::feedSIL(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
1012 std::size_t it0Cnt(0);
1013 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,it0Cnt++)
1015 vtkIdType InfoOnTSId(sil->AddChild(root,edge));
1016 names.push_back((*it0)[0][0].getHumanReadableOverviewOfTS());
1018 vtkIdType NbOfTSId(sil->AddChild(InfoOnTSId,edge));
1019 std::vector<double> ts;
1020 std::vector< std::pair<int,int> > dtits((*it0)[0][0].getTimeStepsInCoarseMEDFileFormat(ts));
1021 std::size_t nbOfTS(dtits.size());
1022 std::ostringstream oss3; oss3 << nbOfTS;
1023 names.push_back(oss3.str());
1024 for(std::size_t i=0;i<nbOfTS;i++)
1026 std::ostringstream oss4; oss4 << dtits[i].first;
1027 vtkIdType DtId(sil->AddChild(NbOfTSId,edge));
1028 names.push_back(oss4.str());
1029 std::ostringstream oss5; oss5 << dtits[i].second;
1030 vtkIdType ItId(sil->AddChild(DtId,edge));
1031 names.push_back(oss5.str());
1032 std::ostringstream oss6; oss6 << ts[i];
1033 sil->AddChild(ItId,edge);
1034 names.push_back(oss6.str());
1037 std::ostringstream oss; oss << MEDFileFieldRepresentationLeavesArrays::TS_STR << it0Cnt;
1038 std::string tsName(oss.str());
1039 vtkIdType typeId0(sil->AddChild(root,edge));
1040 names.push_back(tsName);
1041 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1043 std::string meshName((*it1)[0].getMeshName());
1044 vtkIdType typeId1(sil->AddChild(typeId0,edge));
1045 names.push_back(meshName);
1046 std::size_t it2Cnt(0);
1047 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,it2Cnt++)
1049 std::ostringstream oss2; oss2 << MEDFileFieldRepresentationLeavesArrays::COM_SUP_STR << it2Cnt;
1050 std::string comSupStr(oss2.str());
1051 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1052 names.push_back(comSupStr);
1053 (*it2).feedSIL(_ms,meshName,sil,typeId2,edge,names);
1059 std::string MEDFileFieldRepresentationTree::getActiveMeshName() const
1061 int dummy0(0),dummy1(0),dummy2(0);
1062 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
1063 return leaf.getMeshName();
1066 std::string MEDFileFieldRepresentationTree::feedSILForFamsAndGrps(vtkMutableDirectedGraph* sil, vtkIdType root, vtkVariantArray *edge, std::vector<std::string>& names) const
1068 int dummy0(0),dummy1(0),dummy2(0);
1069 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(dummy0,dummy1,dummy2));
1070 std::string ret(leaf.getMeshName());
1073 for(;i<_ms->getNumberOfMeshes();i++)
1075 m=_ms->getMeshAtPos(i);
1076 if(m->getName()==ret)
1079 if(i==_ms->getNumberOfMeshes())
1080 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::feedSILForFamsAndGrps : internal error #0 !");
1081 vtkIdType typeId0(sil->AddChild(root,edge));
1082 names.push_back(m->getName());
1084 vtkIdType typeId1(sil->AddChild(typeId0,edge));
1085 names.push_back(std::string(ROOT_OF_GRPS_IN_TREE));
1086 std::vector<std::string> grps(m->getGroupsNames());
1087 for(std::vector<std::string>::const_iterator it0=grps.begin();it0!=grps.end();it0++)
1089 vtkIdType typeId2(sil->AddChild(typeId1,edge));
1090 names.push_back(*it0);
1091 std::vector<std::string> famsOnGrp(m->getFamiliesOnGroup((*it0).c_str()));
1092 for(std::vector<std::string>::const_iterator it1=famsOnGrp.begin();it1!=famsOnGrp.end();it1++)
1094 sil->AddChild(typeId2,edge);
1095 names.push_back((*it1).c_str());
1099 vtkIdType typeId11(sil->AddChild(typeId0,edge));
1100 names.push_back(std::string(ROOT_OF_FAM_IDS_IN_TREE));
1101 std::vector<std::string> fams(m->getFamiliesNames());
1102 for(std::vector<std::string>::const_iterator it00=fams.begin();it00!=fams.end();it00++)
1104 sil->AddChild(typeId11,edge);
1105 int famId(m->getFamilyId((*it00).c_str()));
1106 std::ostringstream oss; oss << (*it00) << MEDFileFieldRepresentationLeavesArrays::ZE_SEP << famId;
1107 names.push_back(oss.str());
1112 const MEDFileFieldRepresentationLeavesArrays& MEDFileFieldRepresentationTree::getLeafArr(int id) const
1114 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1115 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1116 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1117 if((*it2).containId(id))
1118 return (*it2).getLeafArr(id);
1119 throw INTERP_KERNEL::Exception("Internal error in MEDFileFieldRepresentationTree::getLeafArr !");
1122 std::string MEDFileFieldRepresentationTree::getNameOf(int id) const
1124 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1125 return elt.getZeName();
1128 const char *MEDFileFieldRepresentationTree::getNameOfC(int id) const
1130 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1131 return elt.getZeNameC();
1134 bool MEDFileFieldRepresentationTree::getStatusOf(int id) const
1136 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1137 return elt.getStatus();
1140 int MEDFileFieldRepresentationTree::getIdHavingZeName(const char *name) const
1143 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1144 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1145 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1146 if((*it2).containZeName(name,ret))
1148 std::ostringstream msg; msg << "MEDFileFieldRepresentationTree::getIdHavingZeName : No such a name \"" << name << "\" !";
1149 throw INTERP_KERNEL::Exception(msg.str().c_str());
1152 bool MEDFileFieldRepresentationTree::changeStatusOfAndUpdateToHaveCoherentVTKDataSet(int id, bool status) const
1154 const MEDFileFieldRepresentationLeavesArrays& elt(getLeafArr(id));
1155 bool ret(elt.setStatus(status));//to be implemented
1159 int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
1162 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1163 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1164 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1165 ret=std::max(ret,(*it2).getNumberOfTS());
1172 void MEDFileFieldRepresentationTree::loadInMemory(MEDCoupling::MEDFileFields *fields, MEDCoupling::MEDFileMeshes *meshes)
1174 _fields=fields; _ms=meshes;
1175 if(_fields.isNotNull())
1180 if(_fields.isNull())
1182 _fields=BuildFieldFromMeshes(_ms);
1186 AppendFieldFromMeshes(_ms,_fields);
1188 _ms->cartesianizeMe();
1189 _fields->removeFieldsWithoutAnyTimeStep();
1190 std::vector<std::string> meshNames(_ms->getMeshesNames());
1191 std::vector< MCAuto<MEDFileFields> > fields_per_mesh(meshNames.size());
1192 for(std::size_t i=0;i<meshNames.size();i++)
1194 fields_per_mesh[i]=_fields->partOfThisLyingOnSpecifiedMeshName(meshNames[i].c_str());
1196 std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS > > allFMTSLeavesToDisplaySafe;
1198 for(std::vector< MCAuto<MEDFileFields> >::const_iterator fields=fields_per_mesh.begin();fields!=fields_per_mesh.end();fields++)
1200 for(int j=0;j<(*fields)->getNumberOfFields();j++)
1202 MCAuto<MEDFileAnyTypeFieldMultiTS> fmts((*fields)->getFieldAtPos((int)j));
1203 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > tmp(fmts->splitDiscretizations());
1205 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::const_iterator it=tmp.begin();it!=tmp.end();it++)
1207 if(!(*it)->presenceOfMultiDiscPerGeoType())
1208 allFMTSLeavesToDisplaySafe.push_back(*it);
1210 {// The case of some parts of field have more than one discretization per geo type.
1211 std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > > subTmp((*it)->splitMultiDiscrPerGeoTypes());
1212 std::size_t it0Cnt(0);
1213 for(std::vector< MCAuto< MEDFileAnyTypeFieldMultiTS > >::iterator it0=subTmp.begin();it0!=subTmp.end();it0++,it0Cnt++)//not const because setName
1215 std::ostringstream oss; oss << (*it0)->getName() << "_" << std::setfill('M') << std::setw(3) << it0Cnt;
1216 (*it0)->setName(oss.str());
1217 allFMTSLeavesToDisplaySafe.push_back(*it0);
1224 std::vector< MEDFileAnyTypeFieldMultiTS *> allFMTSLeavesToDisplay(allFMTSLeavesToDisplaySafe.size());
1225 for(std::size_t i=0;i<allFMTSLeavesToDisplaySafe.size();i++)
1227 allFMTSLeavesToDisplay[i]=allFMTSLeavesToDisplaySafe[i];
1229 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > allFMTSLeavesPerTimeSeries(MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries(allFMTSLeavesToDisplay));
1230 // memory safety part
1231 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > allFMTSLeavesPerTimeSeriesSafe(allFMTSLeavesPerTimeSeries.size());
1232 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeries.size();j++)
1234 allFMTSLeavesPerTimeSeriesSafe[j].resize(allFMTSLeavesPerTimeSeries[j].size());
1235 for(std::size_t k=0;k<allFMTSLeavesPerTimeSeries[j].size();k++)
1237 allFMTSLeavesPerTimeSeries[j][k]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitIntoCommonTimeSeries do not increments the counter
1238 allFMTSLeavesPerTimeSeriesSafe[j][k]=allFMTSLeavesPerTimeSeries[j][k];
1241 // end of memory safety part
1242 // 1st : timesteps, 2nd : meshName, 3rd : common support
1243 this->_data_structure.resize(allFMTSLeavesPerTimeSeriesSafe.size());
1244 for(std::size_t i=0;i<allFMTSLeavesPerTimeSeriesSafe.size();i++)
1246 std::vector< std::string > meshNamesLoc;
1247 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > splitByMeshName;
1248 for(std::size_t j=0;j<allFMTSLeavesPerTimeSeriesSafe[i].size();j++)
1250 std::string meshName(allFMTSLeavesPerTimeSeriesSafe[i][j]->getMeshName());
1251 std::vector< std::string >::iterator it(std::find(meshNamesLoc.begin(),meshNamesLoc.end(),meshName));
1252 if(it==meshNamesLoc.end())
1254 meshNamesLoc.push_back(meshName);
1255 splitByMeshName.resize(splitByMeshName.size()+1);
1256 splitByMeshName.back().push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1259 splitByMeshName[std::distance(meshNamesLoc.begin(),it)].push_back(allFMTSLeavesPerTimeSeriesSafe[i][j]);
1261 _data_structure[i].resize(meshNamesLoc.size());
1262 for(std::size_t j=0;j<splitByMeshName.size();j++)
1264 std::vector< MCAuto<MEDFileFastCellSupportComparator> > fsp;
1265 std::vector< MEDFileAnyTypeFieldMultiTS *> sbmn(splitByMeshName[j].size());
1266 for(std::size_t k=0;k<splitByMeshName[j].size();k++)
1267 sbmn[k]=splitByMeshName[j][k];
1268 //getMeshWithName does not return a newly allocated object ! It is a true get* method !
1269 std::vector< std::vector<MEDFileAnyTypeFieldMultiTS *> > commonSupSplit(MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport(sbmn,_ms->getMeshWithName(meshNamesLoc[j].c_str()),fsp));
1270 std::vector< std::vector< MCAuto<MEDFileAnyTypeFieldMultiTS> > > commonSupSplitSafe(commonSupSplit.size());
1271 this->_data_structure[i][j].resize(commonSupSplit.size());
1272 for(std::size_t k=0;k<commonSupSplit.size();k++)
1274 commonSupSplitSafe[k].resize(commonSupSplit[k].size());
1275 for(std::size_t l=0;l<commonSupSplit[k].size();l++)
1277 commonSupSplit[k][l]->incrRef();//because MEDFileAnyTypeFieldMultiTS::SplitPerCommonSupport does not increment pointers !
1278 commonSupSplitSafe[k][l]=commonSupSplit[k][l];
1281 for(std::size_t k=0;k<commonSupSplit.size();k++)
1282 this->_data_structure[i][j][k]=MEDFileFieldRepresentationLeaves(commonSupSplitSafe[k],fsp[k]);
1285 this->removeEmptyLeaves();
1287 this->computeFullNameInLeaves();
1290 void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts)
1292 MCAuto<MEDFileMeshes> ms;
1293 MCAuto<MEDFileFields> fields;
1296 if((iPart==-1 && nbOfParts==-1) || (iPart==0 && nbOfParts==1))
1298 MCAuto<MEDFileMeshSupports> msups(MEDFileMeshSupports::New(fileName));
1299 MCAuto<MEDFileStructureElements> mse(MEDFileStructureElements::New(fileName,msups));
1300 ms=MEDFileMeshes::New(fileName);
1301 fields=MEDFileFields::NewWithDynGT(fileName,mse,false);//false is important to not read the values
1302 if(ms->presenceOfStructureElements())
1304 fields->loadArrays();
1305 fields->blowUpSE(ms,mse);
1307 int nbMeshes(ms->getNumberOfMeshes());
1308 for(int i=0;i<nbMeshes;i++)
1310 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1311 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1313 tmp2->forceComputationOfParts();
1318 #ifdef MEDREADER_USE_MPI
1319 ms=ParaMEDFileMeshes::New(iPart,nbOfParts,fileName);
1320 int nbMeshes(ms->getNumberOfMeshes());
1321 for(int i=0;i<nbMeshes;i++)
1323 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1324 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1326 MCAuto<DataArrayIdType> tmp3(tmp2->zipCoords());
1328 fields=MEDFileFields::LoadPartOf(fileName,false,ms);//false is important to not read the values
1330 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::loadMainStructureOfFile : request for iPart/nbOfParts=" << iPart << "/" << nbOfParts << " whereas Plugin not compiled with MPI !";
1331 throw INTERP_KERNEL::Exception(oss.str().c_str());
1337 MCAuto<MEDCoupling::SauvReader> sr(MEDCoupling::SauvReader::New(fileName));
1338 MCAuto<MEDCoupling::MEDFileData> mfd(sr->loadInMEDFileDS());
1339 ms=mfd->getMeshes(); ms->incrRef();
1340 int nbMeshes(ms->getNumberOfMeshes());
1341 for(int i=0;i<nbMeshes;i++)
1343 MEDCoupling::MEDFileMesh *tmp(ms->getMeshAtPos(i));
1344 MEDCoupling::MEDFileUMesh *tmp2(dynamic_cast<MEDCoupling::MEDFileUMesh *>(tmp));
1346 tmp2->forceComputationOfParts();
1348 fields=mfd->getFields();
1349 if(fields.isNotNull())
1352 loadInMemory(fields,ms);
1355 void MEDFileFieldRepresentationTree::removeEmptyLeaves()
1357 std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > > newSD;
1358 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1360 std::vector< std::vector< MEDFileFieldRepresentationLeaves > > newSD0;
1361 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1363 std::vector< MEDFileFieldRepresentationLeaves > newSD1;
1364 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1366 newSD1.push_back(*it2);
1368 newSD0.push_back(newSD1);
1371 newSD.push_back(newSD0);
1375 bool MEDFileFieldRepresentationTree::IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos)
1377 if(compInfos.size()!=1)
1379 return compInfos[0]==COMPO_STR_TO_LOCATE_MESH_DA;
1382 std::string MEDFileFieldRepresentationTree::getDftMeshName() const
1384 return _data_structure[0][0][0].getMeshName();
1387 std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, const TimeKeeper& tk) const
1390 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1391 return leaf.getTimeSteps(tk);
1394 vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, ExportedTinyInfo *internalInfo) const
1397 const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
1398 meshName=leaf.getMeshName();
1399 std::vector<double> ts(leaf.getTimeSteps(tk));
1400 std::size_t zeTimeId(0);
1403 std::vector<double> ts2(ts.size());
1404 std::transform(ts.begin(),ts.end(),ts2.begin(),std::bind2nd(std::plus<double>(),-timeReq));
1405 std::transform(ts2.begin(),ts2.end(),ts2.begin(),std::ptr_fun<double,double>(fabs));
1406 zeTimeId=std::distance(ts2.begin(),std::find_if(ts2.begin(),ts2.end(),std::bind2nd(std::less<double>(),1e-14)));
1409 if(zeTimeId==ts.size())
1410 zeTimeId=std::distance(ts.begin(),std::find(ts.begin(),ts.end(),timeReq));
1411 if(zeTimeId==ts.size())
1412 {//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.
1413 //In this case the default behaviour is taken. Keep the highest time step in this lower than timeReq.
1415 double valAttachedToPos(-std::numeric_limits<double>::max());
1416 for(std::size_t i=0;i<ts.size();i++)
1420 if(ts[i]>valAttachedToPos)
1423 valAttachedToPos=ts[i];
1428 {// timeReq is lower than all time steps (ts). So let's keep the lowest time step greater than timeReq.
1429 valAttachedToPos=std::numeric_limits<double>::max();
1430 for(std::size_t i=0;i<ts.size();i++)
1432 if(ts[i]<valAttachedToPos)
1435 valAttachedToPos=ts[i];
1440 std::ostringstream oss; oss.precision(15); oss << "request for time " << timeReq << " but not in ";
1441 std::copy(ts.begin(),ts.end(),std::ostream_iterator<double>(oss,","));
1442 oss << " ! Keep time " << valAttachedToPos << " at pos #" << zeTimeId;
1443 std::cerr << oss.str() << std::endl;
1447 tr=new MEDStdTimeReq((int)zeTimeId);
1449 tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
1450 vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms,internalInfo));
1455 const MEDFileFieldRepresentationLeaves& MEDFileFieldRepresentationTree::getTheSingleActivated(int& lev0, int& lev1, int& lev2) const
1457 int nbOfActivated(0);
1458 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1459 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1460 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1461 if((*it2).isActivated())
1463 if(nbOfActivated!=1)
1465 std::ostringstream oss; oss << "MEDFileFieldRepresentationTree::getTheSingleActivated : Only one leaf must be activated ! Having " << nbOfActivated << " !";
1466 throw INTERP_KERNEL::Exception(oss.str().c_str());
1468 int i0(0),i1(0),i2(0);
1469 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i0++)
1470 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,i1++)
1471 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,i2++)
1472 if((*it2).isActivated())
1474 lev0=i0; lev1=i1; lev2=i2;
1477 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::getTheSingleActivated : Internal error !");
1480 void MEDFileFieldRepresentationTree::printMySelf(std::ostream& os) const
1483 os << "#############################################" << std::endl;
1484 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++,i++)
1487 os << "TS" << i << std::endl;
1488 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++,j++)
1491 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++,k++)
1494 os << " " << (*it2).getMeshName() << std::endl;
1495 os << " Comp" << k << std::endl;
1496 (*it2).printMySelf(os);
1500 os << "$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$" << std::endl;
1503 std::map<std::string,bool> MEDFileFieldRepresentationTree::dumpState() const
1505 std::map<std::string,bool> ret;
1506 for(std::vector< std::vector< std::vector< MEDFileFieldRepresentationLeaves > > >::const_iterator it0=_data_structure.begin();it0!=_data_structure.end();it0++)
1507 for(std::vector< std::vector< MEDFileFieldRepresentationLeaves > >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
1508 for(std::vector< MEDFileFieldRepresentationLeaves >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
1509 (*it2).dumpState(ret);
1513 void MEDFileFieldRepresentationTree::AppendFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms, MEDCoupling::MEDFileFields *ret)
1516 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::AppendFieldFromMeshes : internal error ! NULL ret !");
1517 for(int i=0;i<ms->getNumberOfMeshes();i++)
1519 MEDFileMesh *mm(ms->getMeshAtPos(i));
1520 std::vector<int> levs(mm->getNonEmptyLevels());
1521 MEDCoupling::MCAuto<MEDCoupling::MEDFileField1TS> f1tsMultiLev(MEDCoupling::MEDFileField1TS::New());
1522 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
1525 for(std::vector<int>::const_iterator it=levs.begin();it!=levs.end();it++)
1527 std::vector<INTERP_KERNEL::NormalizedCellType> gts(mmu->getGeoTypesAtLevel(*it));
1528 for(std::vector<INTERP_KERNEL::NormalizedCellType>::const_iterator gt=gts.begin();gt!=gts.end();gt++)
1530 MEDCoupling::MEDCouplingMesh *m(mmu->getDirectUndergroundSingleGeoTypeMesh(*gt));
1531 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1533 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1534 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1537 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1538 f1tsMultiLev->setFieldNoProfileSBT(f);
1543 std::vector<int> levsExt(mm->getNonEmptyLevelsExt());
1544 if(levsExt.size()==levs.size()+1)
1546 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(1));
1547 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_NODES));
1549 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(m->getNumberOfNodes());
1550 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1551 arr->iota(); f->setArray(arr);
1552 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1553 f1tsMultiLev->setFieldNoProfileSBT(f);
1561 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> m(mm->getMeshAtLevel(0));
1562 MEDCoupling::MCAuto<MEDCoupling::MEDCouplingFieldDouble> f(MEDCoupling::MEDCouplingFieldDouble::New(MEDCoupling::ON_CELLS));
1564 MEDCoupling::MCAuto<MEDCoupling::DataArrayDouble> arr(MEDCoupling::DataArrayDouble::New()); arr->alloc(f->getNumberOfTuplesExpected());
1565 arr->setInfoOnComponent(0,std::string(COMPO_STR_TO_LOCATE_MESH_DA));
1568 f->setName(BuildAUniqueArrayNameForMesh(mm->getName(),ret));
1569 f1tsMultiLev->setFieldNoProfileSBT(f);
1572 MEDCoupling::MCAuto<MEDCoupling::MEDFileFieldMultiTS> fmtsMultiLev(MEDCoupling::MEDFileFieldMultiTS::New());
1573 fmtsMultiLev->pushBackTimeStep(f1tsMultiLev);
1574 ret->pushField(fmtsMultiLev);
1578 std::string MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh(const std::string& meshName, const MEDCoupling::MEDFileFields *ret)
1580 const char KEY_STR_TO_AVOID_COLLIDE[]="MESH@";
1582 throw INTERP_KERNEL::Exception("MEDFileFieldRepresentationTree::BuildAUniqueArrayNameForMesh : internal error ! NULL ret !");
1583 std::vector<std::string> fieldNamesAlreadyExisting(ret->getFieldsNames());
1584 if(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),meshName)==fieldNamesAlreadyExisting.end())
1586 std::string tmpName(KEY_STR_TO_AVOID_COLLIDE); tmpName+=meshName;
1587 while(std::find(fieldNamesAlreadyExisting.begin(),fieldNamesAlreadyExisting.end(),tmpName)!=fieldNamesAlreadyExisting.end())
1588 tmpName=std::string(KEY_STR_TO_AVOID_COLLIDE)+tmpName;
1592 MEDCoupling::MEDFileFields *MEDFileFieldRepresentationTree::BuildFieldFromMeshes(const MEDCoupling::MEDFileMeshes *ms)
1594 MEDCoupling::MCAuto<MEDCoupling::MEDFileFields> ret(MEDCoupling::MEDFileFields::New());
1595 AppendFieldFromMeshes(ms,ret);
1599 std::vector<std::string> MEDFileFieldRepresentationTree::SplitFieldNameIntoParts(const std::string& fullFieldName, char sep)
1601 std::vector<std::string> ret;
1603 while(pos!=std::string::npos)
1605 std::size_t curPos(fullFieldName.find_first_of(sep,pos));
1606 std::string elt(fullFieldName.substr(pos,curPos!=std::string::npos?curPos-pos:std::string::npos));
1608 pos=fullFieldName.find_first_not_of(sep,curPos);
1614 * Here the non regression tests.
1615 * const char inp0[]="";
1616 * const char exp0[]="";
1617 * const char inp1[]="field";
1618 * const char exp1[]="field";
1619 * const char inp2[]="_________";
1620 * const char exp2[]="_________";
1621 * const char inp3[]="field_p";
1622 * const char exp3[]="field_p";
1623 * const char inp4[]="field__p";
1624 * const char exp4[]="field_p";
1625 * const char inp5[]="field_p__";
1626 * const char exp5[]="field_p";
1627 * const char inp6[]="field_p_";
1628 * const char exp6[]="field_p";
1629 * const char inp7[]="field_____EDFGEG//sdkjf_____PP_______________";
1630 * const char exp7[]="field_EDFGEG//sdkjf_PP";
1631 * const char inp8[]="field_____EDFGEG//sdkjf_____PP";
1632 * const char exp8[]="field_EDFGEG//sdkjf_PP";
1633 * const char inp9[]="_field_____EDFGEG//sdkjf_____PP_______________";
1634 * const char exp9[]="field_EDFGEG//sdkjf_PP";
1635 * const char inp10[]="___field_____EDFGEG//sdkjf_____PP_______________";
1636 * const char exp10[]="field_EDFGEG//sdkjf_PP";
1638 std::string MEDFileFieldRepresentationTree::PostProcessFieldName(const std::string& fullFieldName)
1640 const char SEP('_');
1641 std::vector<std::string> v(SplitFieldNameIntoParts(fullFieldName,SEP));
1643 return fullFieldName;//should never happen
1647 return fullFieldName;
1651 std::string ret(v[0]);
1652 for(std::size_t i=1;i<v.size();i++)
1657 { ret+=SEP; ret+=v[i]; }
1663 return fullFieldName;
1669 TimeKeeper::TimeKeeper(int policy):_policy(policy)
1673 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1678 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1680 return getTimeStepsRegardingPolicy0(tsPairs,ts);
1682 throw INTERP_KERNEL::Exception("TimeKeeper::getTimeStepsRegardingPolicy : only policy 0 and 1 supported presently !");
1688 * if all of ts are in -1e299,1e299 and different each other pairs are ignored ts taken directly.
1689 * if all of ts are in -1e299,1e299 but some are not different each other ts are ignored pairs used
1690 * if some of ts are out of -1e299,1e299 ts are ignored pairs used
1692 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy0(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1694 std::size_t sz(ts.size());
1695 bool isInHumanRange(true);
1697 for(std::size_t i=0;i<sz;i++)
1700 if(ts[i]<=-1e299 || ts[i]>=1e299)
1701 isInHumanRange=false;
1704 return processedUsingPairOfIds(tsPairs);
1706 return processedUsingPairOfIds(tsPairs);
1707 _postprocessed_time=ts;
1708 return getPostProcessedTime();
1713 * idem than 0, except that ts is preaccumulated before invoking policy 0.
1715 std::vector<double> TimeKeeper::getTimeStepsRegardingPolicy1(const std::vector< std::pair<int,int> >& tsPairs, const std::vector<double>& ts) const
1717 std::size_t sz(ts.size());
1718 std::vector<double> ts2(sz);
1720 for(std::size_t i=0;i<sz;i++)
1725 return getTimeStepsRegardingPolicy0(tsPairs,ts2);
1728 int TimeKeeper::getTimeStepIdFrom(double timeReq) const
1730 std::size_t pos(std::distance(_postprocessed_time.begin(),std::find(_postprocessed_time.begin(),_postprocessed_time.end(),timeReq)));
1734 void TimeKeeper::printSelf(std::ostream& oss) const
1736 std::size_t sz(_activated_ts.size());
1737 for(std::size_t i=0;i<sz;i++)
1739 oss << "(" << i << "," << _activated_ts[i].first << "), ";
1743 std::vector<bool> TimeKeeper::getTheVectOfBool() const
1745 std::size_t sz(_activated_ts.size());
1746 std::vector<bool> ret(sz);
1747 for(std::size_t i=0;i<sz;i++)
1749 ret[i]=_activated_ts[i].first;
1754 std::vector<double> TimeKeeper::processedUsingPairOfIds(const std::vector< std::pair<int,int> >& tsPairs) const
1756 std::size_t sz(tsPairs.size());
1757 std::set<int> s0,s1;
1758 for(std::size_t i=0;i<sz;i++)
1759 { s0.insert(tsPairs[i].first); s1.insert(tsPairs[i].second); }
1762 _postprocessed_time.resize(sz);
1763 for(std::size_t i=0;i<sz;i++)
1764 _postprocessed_time[i]=(double)tsPairs[i].first;
1765 return getPostProcessedTime();
1769 _postprocessed_time.resize(sz);
1770 for(std::size_t i=0;i<sz;i++)
1771 _postprocessed_time[i]=(double)tsPairs[i].second;
1772 return getPostProcessedTime();
1774 //TimeKeeper::processedUsingPairOfIds : you are not a lucky guy ! All your time steps info in MEDFile are not discriminant taken one by one !
1775 _postprocessed_time.resize(sz);
1776 for(std::size_t i=0;i<sz;i++)
1777 _postprocessed_time[i]=(double)i;
1778 return getPostProcessedTime();
1781 void TimeKeeper::setMaxNumberOfTimeSteps(int maxNumberOfTS)
1783 _activated_ts.resize(maxNumberOfTS);
1784 for(int i=0;i<maxNumberOfTS;i++)
1786 std::ostringstream oss; oss << "000" << i;
1787 _activated_ts[i]=std::pair<bool,std::string>(true,oss.str());