1 // Copyright (C) 2017-2022 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 (EDF R&D)
21 #include "VTKToMEDMem.h"
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkIntArray.h"
25 #include "vtkLongArray.h"
26 #include "vtkLongLongArray.h"
27 #include "vtkCellData.h"
28 #include "vtkPointData.h"
29 #include "vtkFloatArray.h"
30 #include "vtkCellArray.h"
32 #include "vtkStreamingDemandDrivenPipeline.h"
33 #include "vtkInformationDataObjectMetaDataKey.h"
34 #include "vtkUnstructuredGrid.h"
35 #include "vtkMultiBlockDataSet.h"
36 #include "vtkRectilinearGrid.h"
37 #include "vtkInformationStringKey.h"
38 #include "vtkAlgorithmOutput.h"
39 #include "vtkObjectFactory.h"
40 #include "vtkMutableDirectedGraph.h"
41 #include "vtkMultiBlockDataSet.h"
42 #include "vtkPolyData.h"
43 #include "vtkDataSet.h"
44 #include "vtkInformationVector.h"
45 #include "vtkInformation.h"
46 #include "vtkDataArraySelection.h"
47 #include "vtkTimeStamp.h"
48 #include "vtkInEdgeIterator.h"
49 #include "vtkInformationDataObjectKey.h"
50 #include "vtkExecutive.h"
51 #include "vtkVariantArray.h"
52 #include "vtkStringArray.h"
53 #include "vtkDoubleArray.h"
54 #include "vtkCharArray.h"
55 #include "vtkUnsignedCharArray.h"
56 #include "vtkDataSetAttributes.h"
57 #include "vtkDemandDrivenPipeline.h"
58 #include "vtkDataObjectTreeIterator.h"
59 #include "vtkWarpScalar.h"
66 using VTKToMEDMemDevelopSurface::Grp;
67 using VTKToMEDMemDevelopSurface::Fam;
69 using MEDCoupling::MEDFileData;
70 using MEDCoupling::MEDFileMesh;
71 using MEDCoupling::MEDFileCMesh;
72 using MEDCoupling::MEDFileUMesh;
73 using MEDCoupling::MEDFileFields;
74 using MEDCoupling::MEDFileMeshes;
76 using MEDCoupling::MEDFileInt32Field1TS;
77 using MEDCoupling::MEDFileInt64Field1TS;
78 using MEDCoupling::MEDFileField1TS;
79 using MEDCoupling::MEDFileIntFieldMultiTS;
80 using MEDCoupling::MEDFileFieldMultiTS;
81 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
82 using MEDCoupling::DataArray;
83 using MEDCoupling::DataArrayInt32;
84 using MEDCoupling::DataArrayInt64;
85 using MEDCoupling::DataArrayFloat;
86 using MEDCoupling::DataArrayDouble;
87 using MEDCoupling::MEDCouplingMesh;
88 using MEDCoupling::MEDCouplingUMesh;
89 using MEDCoupling::MEDCouplingCMesh;
90 using MEDCoupling::MEDCouplingFieldDouble;
91 using MEDCoupling::MEDCouplingFieldFloat;
92 using MEDCoupling::MEDCouplingFieldInt;
93 using MEDCoupling::MEDCouplingFieldInt64;
94 using MEDCoupling::MCAuto;
95 using MEDCoupling::Traits;
96 using MEDCoupling::MLFieldTraits;
98 using namespace VTKToMEDMemDevelopSurface;
102 Fam::Fam(const std::string& name)
104 static const char ZE_SEP[]="@@][@@";
105 std::size_t pos(name.find(ZE_SEP));
106 std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
107 std::istringstream iss(name1);
114 #include "VTKMEDTraits.hxx"
116 static std::map<int,int> ComputeMapOfType()
118 std::map<int,int> ret;
119 int nbOfTypesInMC(sizeof(MEDCOUPLING2VTKTYPETRADUCER)/sizeof( decltype(MEDCOUPLING2VTKTYPETRADUCER[0]) ));
120 for(int i=0;i<nbOfTypesInMC;i++)
122 auto vtkId(MEDCOUPLING2VTKTYPETRADUCER[i]);
123 if(vtkId!=MEDCOUPLING2VTKTYPETRADUCER_NONE)
129 static std::string GetMeshNameWithContext(const std::vector<int>& context)
131 static const char DFT_MESH_NAME[]="Mesh";
133 return DFT_MESH_NAME;
134 std::ostringstream oss; oss << DFT_MESH_NAME;
135 for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
140 static DataArrayIdType *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
143 throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
144 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
145 std::size_t nbElts(nbTuples*nbComp);
146 MCAuto<DataArrayIdType> ret(DataArrayIdType::New());
147 ret->alloc(nbTuples,nbComp);
148 for(int i=0;i<nbComp;i++)
150 const char *comp(data->GetComponentName(i));
152 ret->setInfoOnComponent(i,comp);
154 mcIdType *ptOut(ret->getPointer());
155 vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
158 const int *pt(d0->GetPointer(0));
159 std::copy(pt,pt+nbElts,ptOut);
162 vtkLongArray *d1(vtkLongArray::SafeDownCast(data));
165 const long *pt(d1->GetPointer(0));
166 std::copy(pt,pt+nbElts,ptOut);
169 vtkLongLongArray *d1l(vtkLongLongArray::SafeDownCast(data));
172 const long long *pt(d1l->GetPointer(0));
173 std::copy(pt,pt+nbElts,ptOut);
176 vtkIdTypeArray *d3(vtkIdTypeArray::SafeDownCast(data));
179 const vtkIdType *pt(d3->GetPointer(0));
180 std::copy(pt,pt+nbElts,ptOut);
183 vtkUnsignedCharArray *d2(vtkUnsignedCharArray::SafeDownCast(data));
186 const unsigned char *pt(d2->GetPointer(0));
187 std::copy(pt,pt+nbElts,ptOut);
190 std::ostringstream oss;
191 oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
192 throw MZCException(oss.str());
196 static typename Traits<T>::ArrayType *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
199 throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
200 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
201 std::size_t nbElts(nbTuples*nbComp);
202 MCAuto< typename Traits<T>::ArrayType > ret(Traits<T>::ArrayType::New());
203 ret->alloc(nbTuples,nbComp);
204 for(int i=0;i<nbComp;i++)
206 const char *comp(data->GetComponentName(i));
208 ret->setInfoOnComponent(i,comp);
211 if(nbComp>1 && nbComp<=3)
214 tmp[0]=(char)('X'+i); tmp[1]='\0';
215 ret->setInfoOnComponent(i,tmp);
219 T *ptOut(ret->getPointer());
220 typename MEDFileVTKTraits<T>::VtkType *d0(MEDFileVTKTraits<T>::VtkType::SafeDownCast(data));
223 const T *pt(d0->GetPointer(0));
224 for(std::size_t i=0;i<nbElts;i++)
228 std::ostringstream oss;
229 oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << data->GetClassName() << "\" type !";
230 throw MZCException(oss.str());
233 static DataArrayDouble *ConvertVTKArrayToMCArrayDoubleForced(vtkDataArray *data)
236 throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : internal error 0 !");
237 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
240 MCAuto<DataArrayFloat> ret(ConvertVTKArrayToMCArrayDouble<float>(data));
241 MCAuto<DataArrayDouble> ret2(ret->convertToDblArr());
244 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
246 return ConvertVTKArrayToMCArrayDouble<double>(data);
247 throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : unrecognized type of data for double !");
250 static DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
253 throw MZCException("ConvertVTKArrayToMCArray : internal error !");
254 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
256 return ConvertVTKArrayToMCArrayDouble<float>(data);
257 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
259 return ConvertVTKArrayToMCArrayDouble<double>(data);
260 vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
261 vtkLongArray *d3(vtkLongArray::SafeDownCast(data));
262 vtkUnsignedCharArray *d4(vtkUnsignedCharArray::SafeDownCast(data));
263 vtkIdTypeArray *d5(vtkIdTypeArray::SafeDownCast(data));
264 vtkLongLongArray *d6(vtkLongLongArray::SafeDownCast(data));
265 if(d2 || d3 || d4 || d5 || d6)
266 return ConvertVTKArrayToMCArrayInt(data);
267 std::ostringstream oss;
268 oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
269 throw MZCException(oss.str());
272 static MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
274 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
275 subMesh->setCoords(coords); subMesh->allocateCells();
276 int nbCells(ca->GetNumberOfCells());
279 //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
280 const vtkIdType *conn(ca->GetData()->GetPointer(0));
281 for(int i=0;i<nbCells;i++)
283 mcIdType sz(ToIdType(*conn++));
284 std::vector<mcIdType> conn2(sz);
285 for(int jj=0;jj<sz;jj++)
286 conn2[jj]=ToIdType(conn[jj]);
287 subMesh->insertNextCell(type,sz,&conn2[0]);
290 return subMesh.retn();
293 static MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayIdType>& ids)
295 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
296 subMesh->setCoords(coords); subMesh->allocateCells();
297 int nbCells(ca->GetNumberOfCells());
300 //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
301 const vtkIdType *conn(ca->GetData()->GetPointer(0));
302 ids=DataArrayIdType::New() ; ids->alloc(0,1);
303 for(int i=0;i<nbCells;i++)
309 for(int j=0;j<nbTri;j++,conn++)
311 mcIdType conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
312 subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
313 ids->pushBackSilent(i);
318 std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
319 throw MZCException(oss.str());
323 return subMesh.retn();
329 MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
330 MicroField(const std::vector< MicroField >& vs);
331 void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
332 MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
333 std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
335 MCAuto<MEDCouplingUMesh> _m;
336 std::vector<MCAuto<DataArray> > _cellFs;
337 std::vector<MCAuto<DataArray> > _nodeFs;
340 MicroField::MicroField(const std::vector< MicroField >& vs)
342 std::size_t sz(vs.size());
343 std::vector<const MEDCouplingUMesh *> vs2(sz);
344 std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
346 for(std::size_t ii=0;ii<sz;ii++)
348 vs2[ii]=vs[ii].getMesh();
349 arrs2[ii]=vs[ii].getCellFields();
351 nbElts=(int)arrs2[ii].size();
353 if((int)arrs2[ii].size()!=nbElts)
354 throw MZCException("MicroField cstor : internal error !");
356 _cellFs.resize(nbElts);
357 for(int ii=0;ii<nbElts;ii++)
359 std::vector<const DataArray *> arrsTmp(sz);
360 for(std::size_t jj=0;jj<sz;jj++)
362 arrsTmp[jj]=arrs2[jj][ii];
364 _cellFs[ii]=DataArray::Aggregate(arrsTmp);
366 _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
370 static void AppendToFields(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, const DataArrayIdType *n2oPtr, typename MEDFileVTKTraits<T>::MCType *dadPtr, MEDFileFields *fs, double timeStep, int tsId)
372 std::string fieldName(dadPtr->getName());
373 MCAuto< typename Traits<T>::FieldType > f(Traits<T>::FieldType::New(tf));
374 f->setTime(timeStep,tsId,0);
376 std::string fieldNameForChuckNorris(MEDCoupling::MEDFileAnyTypeField1TSWithoutSDA::FieldNameToMEDFileConvention(fieldName));
377 f->setName(fieldNameForChuckNorris);
383 MCAuto< typename Traits<T>::ArrayType > dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
387 MCAuto< typename MLFieldTraits<T>::FMTSType > fmts(MLFieldTraits<T>::FMTSType::New());
388 MCAuto< typename MLFieldTraits<T>::F1TSType > f1ts(MLFieldTraits<T>::F1TSType::New());
389 f1ts->setFieldNoProfileSBT(f);
390 fmts->pushBackTimeStep(f1ts);
394 static void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayIdType *n2oPtr, double timeStep, int tsId)
396 static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
397 static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
398 if(!da || !mesh || !mfd)
399 throw MZCException("AppendMCFieldFrom : internal error !");
400 MEDFileFields *fs(mfd->getFields());
401 MEDFileMeshes *ms(mfd->getMeshes());
403 throw MZCException("AppendMCFieldFrom : internal error 2 !");
404 MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
407 AppendToFields<double>(tf,mesh,n2oPtr,dad,fs,timeStep,tsId);
410 MCAuto<DataArrayFloat> daf(MEDCoupling::DynamicCast<DataArray,DataArrayFloat>(da));
413 AppendToFields<float>(tf,mesh,n2oPtr,daf,fs,timeStep,tsId);
416 MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
417 MCAuto<DataArrayIdType> daId(MEDCoupling::DynamicCast<DataArray,DataArrayIdType>(da));
418 if(dai.isNotNull() || daId.isNotNull())
420 std::string fieldName(da->getName());
421 if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
424 AppendToFields<mcIdType>(tf,mesh,n2oPtr,daId,fs,timeStep,tsId);
426 AppendToFields<int>(tf,mesh,n2oPtr,dai,fs,timeStep,tsId);
429 else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
431 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
433 throw MZCException("AppendMCFieldFrom : internal error 3 !");
435 throw MZCException("AppendMCFieldFrom : internal error 3 (not mcIdType) !");
437 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daId);
440 MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
441 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
444 else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
446 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
448 throw MZCException("AppendMCFieldFrom : internal error 4 !");
450 throw MZCException("AppendMCFieldFrom : internal error 4 (not mcIdType) !");
452 mm->setFamilyFieldArr(1,daId);
455 MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
456 mm->setFamilyFieldArr(1,dai2);
462 static void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf, double timeStep, int tsId)
465 throw MZCException("PutAtLevelDealOrder : internal error !");
466 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
467 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
469 throw MZCException("PutAtLevelDealOrder : internal error 2 !");
470 MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
471 mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
472 MCAuto<DataArrayIdType> o2n(mesh->sortCellsInMEDFileFrmt());
473 //const DataArrayIdType *o2nPtr(o2n); // todo: unused
474 MCAuto<DataArrayIdType> n2o;
475 mmu->setMeshAtLevel(meshDimRel,mesh);
476 const DataArrayIdType *n2oPtr(0);
479 n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
481 if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
484 mm->setRenumFieldArr(meshDimRel,n2o);
487 std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
488 for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
490 MCAuto<DataArray> da(*it);
491 AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr,timeStep,tsId);
495 static void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms, double timeStep, int tsId)
498 throw MZCException("AssignSingleGTMeshes : internal error !");
499 MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
500 MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
502 throw MZCException("AssignSingleGTMeshes : internal error 2 !");
503 int meshDim(-std::numeric_limits<int>::max());
504 std::map<int, std::vector< MicroField > > ms2;
505 for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
507 const MEDCouplingUMesh *elt((*it).getMesh());
510 int myMeshDim(elt->getMeshDimension());
511 meshDim=std::max(meshDim,myMeshDim);
512 ms2[myMeshDim].push_back(*it);
517 for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
519 const std::vector< MicroField >& vs((*it).second);
522 PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0],timeStep,tsId);
526 MicroField merge(vs);
527 PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge,timeStep,tsId);
532 static DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
535 throw MZCException("BuildCoordsFrom : internal error !");
536 vtkPoints *pts(ds->GetPoints());
538 throw MZCException("BuildCoordsFrom : internal error 2 !");
539 vtkDataArray *data(pts->GetData());
541 throw MZCException("BuildCoordsFrom : internal error 3 !");
542 return ConvertVTKArrayToMCArrayDoubleForced(data);
545 static void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa, double timeStep, int tsId)
548 throw MZCException("AddNodeFields : internal error !");
549 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
550 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
552 throw MZCException("AddNodeFields : internal error 2 !");
553 MCAuto<MEDCouplingUMesh> mesh;
554 if(!mmu->getNonEmptyLevels().empty())
555 mesh=mmu->getMeshAtLevel(0);
558 mesh=MEDCouplingUMesh::Build0DMeshFromCoords(mmu->getCoords());
559 mesh->setName(mmu->getName());
561 int nba(dsa->GetNumberOfArrays());
562 for(int i=0;i<nba;i++)
564 vtkDataArray *arr(dsa->GetArray(i));
565 const char *name(arr->GetName());
568 MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
570 AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,NULL,timeStep,tsId);
574 static std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayIdType *part, vtkDataSetAttributes *dsa)
576 std::vector< MCAuto<DataArray> > ret;
579 int nba(dsa->GetNumberOfArrays());
580 for(int i=0;i<nba;i++)
582 vtkDataArray *arr(dsa->GetArray(i));
585 const char *name(arr->GetName());
586 //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
587 //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
588 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
590 mcarr=mcarr->selectByTupleId(part->begin(),part->end());
591 mcarr->setName(name);
592 ret.push_back(mcarr);
597 static std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
599 std::vector< MCAuto<DataArray> > ret;
602 int nba(dsa->GetNumberOfArrays());
603 for(int i=0;i<nba;i++)
605 vtkDataArray *arr(dsa->GetArray(i));
608 const char *name(arr->GetName());
609 //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
610 //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
611 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
612 mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
613 mcarr->setName(name);
614 ret.push_back(mcarr);
619 static void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
622 throw MZCException("ConvertFromRectilinearGrid : internal error !");
624 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
625 ret->setMeshes(meshes);
626 MCAuto<MEDFileFields> fields(MEDFileFields::New());
627 ret->setFields(fields);
629 MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
630 meshes->pushMesh(cmesh);
631 MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
632 vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
635 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cx));
636 cmeshmc->setCoordsAt(0,arr);
640 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cy));
641 cmeshmc->setCoordsAt(1,arr);
645 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cz));
646 cmeshmc->setCoordsAt(2,arr);
648 std::string meshName(GetMeshNameWithContext(context));
649 cmeshmc->setName(meshName);
650 cmesh->setMesh(cmeshmc);
651 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
652 for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
654 MCAuto<DataArray> da(*it);
655 AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,NULL,timeStep,tsId);
657 std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
658 for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
660 MCAuto<DataArray> da(*it);
661 AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,NULL,timeStep,tsId);
665 static void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context, double timeStep, int tsId)
668 throw MZCException("ConvertFromPolyData : internal error !");
670 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
671 ret->setMeshes(meshes);
672 MCAuto<MEDFileFields> fields(MEDFileFields::New());
673 ret->setFields(fields);
675 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
676 meshes->pushMesh(umesh);
677 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
678 umesh->setCoords(coords);
679 umesh->setName(GetMeshNameWithContext(context));
682 std::vector< MicroField > ms;
683 vtkCellArray *cd(ds->GetVerts());
686 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
687 if((const MEDCouplingUMesh *)subMesh)
689 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
690 offset+=subMesh->getNumberOfCells();
691 ms.push_back(MicroField(subMesh,cellFs));
694 vtkCellArray *cc(ds->GetLines());
697 MCAuto<MEDCouplingUMesh> subMesh;
700 subMesh=BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2);
702 catch(INTERP_KERNEL::Exception& e)
704 std::ostringstream oss; oss << "MEDWriter does not manage polyline cell type because MED file format does not support it ! Maybe it is the source of the problem ? The cause of this exception was " << e.what() << std::endl;
705 throw INTERP_KERNEL::Exception(oss.str());
707 if((const MEDCouplingUMesh *)subMesh)
709 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
710 offset+=subMesh->getNumberOfCells();
711 ms.push_back(MicroField(subMesh,cellFs));
714 vtkCellArray *cb(ds->GetPolys());
717 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
718 if((const MEDCouplingUMesh *)subMesh)
720 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
721 offset+=subMesh->getNumberOfCells();
722 ms.push_back(MicroField(subMesh,cellFs));
725 vtkCellArray *ca(ds->GetStrips());
728 MCAuto<DataArrayIdType> ids;
729 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
730 if((const MEDCouplingUMesh *)subMesh)
732 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
733 offset+=subMesh->getNumberOfCells();
734 ms.push_back(MicroField(subMesh,cellFs));
737 AssignSingleGTMeshes(ret,ms,timeStep,tsId);
738 AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
741 static void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
744 throw MZCException("ConvertFromUnstructuredGrid : internal error !");
746 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
747 ret->setMeshes(meshes);
748 MCAuto<MEDFileFields> fields(MEDFileFields::New());
749 ret->setFields(fields);
751 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
752 meshes->pushMesh(umesh);
753 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
754 umesh->setCoords(coords);
755 umesh->setName(GetMeshNameWithContext(context));
756 vtkIdType nbCells(ds->GetNumberOfCells());
757 vtkCellArray *ca(ds->GetCells());
760 //vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries()); // todo: unused
761 //vtkIdType *caPtr(ca->GetData()->GetPointer(0)); // todo: unused
762 vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
764 throw MZCException("ConvertFromUnstructuredGrid : internal error");
765 vtkIdTypeArray *cla(ds->GetCellLocationsArray());
766 //const vtkIdType *claPtr(cla->GetPointer(0)); // todo: unused
768 throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
769 const unsigned char *ctPtr(ct->GetPointer(0));
770 std::map<int,int> m(ComputeMapOfType());
771 MCAuto<DataArrayInt> lev(DataArrayInt::New()) ; lev->alloc(nbCells,1);
772 int *levPtr(lev->getPointer());
773 for(vtkIdType i=0;i<nbCells;i++)
775 std::map<int,int>::iterator it(m.find(ctPtr[i]));
778 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
779 levPtr[i]=cm.getDimension();
783 if(ctPtr[i]==VTK_POLY_VERTEX)
785 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(INTERP_KERNEL::NORM_POINT1));
786 levPtr[i]=cm.getDimension();
790 std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
791 throw MZCException(oss.str());
795 MCAuto<DataArrayInt> levs(lev->getDifferentValues());
796 std::vector< MicroField > ms;
797 //vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations()); // todo: unused
798 for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
800 MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
801 m0->setCoords(coords); m0->allocateCells();
802 MCAuto<DataArrayIdType> cellIdsCurLev(lev->findIdsEqual(*curLev));
803 for(const mcIdType *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
805 int vtkType(ds->GetCellType(*cellId));
806 std::map<int,int>::iterator it(m.find(vtkType));
807 INTERP_KERNEL::NormalizedCellType ct=it!=m.end()?(INTERP_KERNEL::NormalizedCellType)((*it).second):INTERP_KERNEL::NORM_POINT1;
808 if(ct!=INTERP_KERNEL::NORM_POLYHED && vtkType!=VTK_POLY_VERTEX)
811 const vtkIdType *pts(nullptr);
812 ds->GetCellPoints(*cellId, sz, pts);
813 std::vector<mcIdType> conn2(pts,pts+sz);
814 m0->insertNextCell(ct,sz,conn2.data());
816 else if(ct==INTERP_KERNEL::NORM_POLYHED)
818 // # de faces du polyèdre
819 vtkIdType nbOfFaces(0);
820 // connectivé des faces (numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...)
821 const vtkIdType *facPtr(nullptr);
822 ds->GetFaceStream(*cellId, nbOfFaces, facPtr);
823 std::vector<mcIdType> conn;
824 for(vtkIdType k=0;k<nbOfFaces;k++)
826 vtkIdType nbOfNodesInFace(*facPtr++);
827 std::copy(facPtr,facPtr+nbOfNodesInFace,std::back_inserter(conn));
830 facPtr+=nbOfNodesInFace;
832 m0->insertNextCell(ct,ToIdType(conn.size()),&conn[0]);
837 const vtkIdType *pts(nullptr);
838 ds->GetCellPoints(*cellId, sz, pts);
840 throw MZCException("ConvertFromUnstructuredGrid : non single poly vertex not managed by MED !");
841 m0->insertNextCell(ct,1,(const mcIdType*)pts);
844 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
845 ms.push_back(MicroField(m0,cellFs));
847 AssignSingleGTMeshes(ret,ms,timeStep,tsId);
848 AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
853 void VTKToMEDMemDevelopSurface::WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
856 throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
857 vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
858 vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
859 vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
862 ConvertFromPolyData(mfd,ds2,context,timeStep,tsId);
866 ConvertFromUnstructuredGrid(mfd,ds3,context,timeStep,tsId);
870 ConvertFromRectilinearGrid(mfd,ds4,context,timeStep,tsId);
873 throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
876 static void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
879 throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
880 int nbBlocks(ds->GetNumberOfBlocks());
881 if(nbBlocks==1 && context.empty())
883 vtkDataObject *uniqueElt(ds->GetBlock(0));
885 throw MZCException("Unique elt in multiblock is NULL !");
886 vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
889 WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context,timeStep,tsId);
893 for(int i=0;i<nbBlocks;i++)
895 vtkDataObject *elt(ds->GetBlock(i));
896 std::vector<int> context2;
897 context2.push_back(i);
900 std::ostringstream oss; oss << "In context ";
901 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
902 oss << " at pos #" << i << " elt is NULL !";
903 throw MZCException(oss.str());
905 vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
908 WriteMEDFileFromVTKDataSet(mfd,elt1,context,timeStep,tsId);
911 vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
914 WriteMEDFileFromVTKMultiBlock(mfd,elt2,context,timeStep,tsId);
917 std::ostringstream oss; oss << "In context ";
918 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
919 oss << " at pos #" << i << " elt not recognized data type !";
920 throw MZCException(oss.str());
924 void VTKToMEDMemDevelopSurface::WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input, double timeStep, int tsId)
927 throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
928 std::vector<int> context;
929 vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
932 WriteMEDFileFromVTKDataSet(mfd,input1,context,timeStep,tsId);
935 vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
938 WriteMEDFileFromVTKMultiBlock(mfd,input2,context,timeStep,tsId);
941 throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
944 void VTKToMEDMemDevelopSurface::PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
950 MEDFileMeshes *meshes(mfd->getMeshes());
953 if(meshes->getNumberOfMeshes()!=1)
955 MEDFileMesh *mm(meshes->getMeshAtPos(0));
958 mm->setName(meshName);
959 for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
960 mm->setFamilyId((*it).getName(),(*it).getID());
961 for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
962 mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
963 MEDFileFields *fields(mfd->getFields());
966 for(int i=0;i<fields->getNumberOfFields();i++)
968 MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
971 fmts->setMeshName(meshName);