1 // Copyright (C) 2017-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 (EDF R&D)
21 #include "VTKToMEDMem.h"
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkIntArray.h"
25 #include "vtkLongArray.h"
26 #include "vtkCellData.h"
27 #include "vtkPointData.h"
28 #include "vtkFloatArray.h"
29 #include "vtkCellArray.h"
31 #include "vtkStreamingDemandDrivenPipeline.h"
32 #include "vtkInformationDataObjectMetaDataKey.h"
33 #include "vtkUnstructuredGrid.h"
34 #include "vtkMultiBlockDataSet.h"
35 #include "vtkRectilinearGrid.h"
36 #include "vtkInformationStringKey.h"
37 #include "vtkAlgorithmOutput.h"
38 #include "vtkObjectFactory.h"
39 #include "vtkMutableDirectedGraph.h"
40 #include "vtkMultiBlockDataSet.h"
41 #include "vtkPolyData.h"
42 #include "vtkDataSet.h"
43 #include "vtkInformationVector.h"
44 #include "vtkInformation.h"
45 #include "vtkDataArraySelection.h"
46 #include "vtkTimeStamp.h"
47 #include "vtkInEdgeIterator.h"
48 #include "vtkInformationDataObjectKey.h"
49 #include "vtkExecutive.h"
50 #include "vtkVariantArray.h"
51 #include "vtkStringArray.h"
52 #include "vtkDoubleArray.h"
53 #include "vtkCharArray.h"
54 #include "vtkUnsignedCharArray.h"
55 #include "vtkDataSetAttributes.h"
56 #include "vtkDemandDrivenPipeline.h"
57 #include "vtkDataObjectTreeIterator.h"
58 #include "vtkWarpScalar.h"
65 using VTKToMEDMem::Grp;
66 using VTKToMEDMem::Fam;
68 using MEDCoupling::MEDFileData;
69 using MEDCoupling::MEDFileMesh;
70 using MEDCoupling::MEDFileCMesh;
71 using MEDCoupling::MEDFileUMesh;
72 using MEDCoupling::MEDFileFields;
73 using MEDCoupling::MEDFileMeshes;
75 using MEDCoupling::MEDFileInt32Field1TS;
76 using MEDCoupling::MEDFileInt64Field1TS;
77 using MEDCoupling::MEDFileField1TS;
78 using MEDCoupling::MEDFileIntFieldMultiTS;
79 using MEDCoupling::MEDFileFieldMultiTS;
80 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
81 using MEDCoupling::DataArray;
82 using MEDCoupling::DataArrayInt32;
83 using MEDCoupling::DataArrayInt64;
84 using MEDCoupling::DataArrayFloat;
85 using MEDCoupling::DataArrayDouble;
86 using MEDCoupling::MEDCouplingMesh;
87 using MEDCoupling::MEDCouplingUMesh;
88 using MEDCoupling::MEDCouplingCMesh;
89 using MEDCoupling::MEDCouplingFieldDouble;
90 using MEDCoupling::MEDCouplingFieldFloat;
91 using MEDCoupling::MEDCouplingFieldInt;
92 using MEDCoupling::MEDCouplingFieldInt64;
93 using MEDCoupling::MCAuto;
94 using MEDCoupling::Traits;
95 using MEDCoupling::MLFieldTraits;
99 Fam::Fam(const std::string& name)
101 static const char ZE_SEP[]="@@][@@";
102 std::size_t pos(name.find(ZE_SEP));
103 std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
104 std::istringstream iss(name1);
111 #include "VTKMEDTraits.hxx"
113 std::map<int,int> ComputeMapOfType()
115 std::map<int,int> ret;
116 int nbOfTypesInMC(sizeof(MEDCOUPLING2VTKTYPETRADUCER)/sizeof( decltype(MEDCOUPLING2VTKTYPETRADUCER[0]) ));
117 for(int i=0;i<nbOfTypesInMC;i++)
119 auto vtkId(MEDCOUPLING2VTKTYPETRADUCER[i]);
120 if(vtkId!=MEDCOUPLING2VTKTYPETRADUCER_NONE)
126 std::string GetMeshNameWithContext(const std::vector<int>& context)
128 static const char DFT_MESH_NAME[]="Mesh";
130 return DFT_MESH_NAME;
131 std::ostringstream oss; oss << DFT_MESH_NAME;
132 for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
137 DataArrayIdType *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
140 throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
141 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
142 std::size_t nbElts(nbTuples*nbComp);
143 MCAuto<DataArrayIdType> ret(DataArrayIdType::New());
144 ret->alloc(nbTuples,nbComp);
145 for(int i=0;i<nbComp;i++)
147 const char *comp(data->GetComponentName(i));
149 ret->setInfoOnComponent(i,comp);
151 mcIdType *ptOut(ret->getPointer());
152 vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
155 const int *pt(d0->GetPointer(0));
156 std::copy(pt,pt+nbElts,ptOut);
159 vtkLongArray *d1(vtkLongArray::SafeDownCast(data));
162 const long *pt(d1->GetPointer(0));
163 std::copy(pt,pt+nbElts,ptOut);
166 vtkIdTypeArray *d3(vtkIdTypeArray::SafeDownCast(data));
169 const vtkIdType *pt(d3->GetPointer(0));
170 std::copy(pt,pt+nbElts,ptOut);
173 vtkUnsignedCharArray *d2(vtkUnsignedCharArray::SafeDownCast(data));
176 const unsigned char *pt(d2->GetPointer(0));
177 std::copy(pt,pt+nbElts,ptOut);
180 std::ostringstream oss;
181 oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
182 throw MZCException(oss.str());
186 typename Traits<T>::ArrayType *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
189 throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
190 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
191 std::size_t nbElts(nbTuples*nbComp);
192 MCAuto< typename Traits<T>::ArrayType > ret(Traits<T>::ArrayType::New());
193 ret->alloc(nbTuples,nbComp);
194 for(int i=0;i<nbComp;i++)
196 const char *comp(data->GetComponentName(i));
198 ret->setInfoOnComponent(i,comp);
201 if(nbComp>1 && nbComp<=3)
204 tmp[0]=(char)('X'+i); tmp[1]='\0';
205 ret->setInfoOnComponent(i,tmp);
209 T *ptOut(ret->getPointer());
210 typename MEDFileVTKTraits<T>::VtkType *d0(MEDFileVTKTraits<T>::VtkType::SafeDownCast(data));
213 const T *pt(d0->GetPointer(0));
214 for(std::size_t i=0;i<nbElts;i++)
218 std::ostringstream oss;
219 oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << data->GetClassName() << "\" type !";
220 throw MZCException(oss.str());
223 DataArrayDouble *ConvertVTKArrayToMCArrayDoubleForced(vtkDataArray *data)
226 throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : internal error 0 !");
227 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
230 MCAuto<DataArrayFloat> ret(ConvertVTKArrayToMCArrayDouble<float>(data));
231 MCAuto<DataArrayDouble> ret2(ret->convertToDblArr());
234 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
236 return ConvertVTKArrayToMCArrayDouble<double>(data);
237 throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : unrecognized type of data for double !");
240 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
243 throw MZCException("ConvertVTKArrayToMCArray : internal error !");
244 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
246 return ConvertVTKArrayToMCArrayDouble<float>(data);
247 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
249 return ConvertVTKArrayToMCArrayDouble<double>(data);
250 vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
251 vtkLongArray *d3(vtkLongArray::SafeDownCast(data));
252 vtkUnsignedCharArray *d4(vtkUnsignedCharArray::SafeDownCast(data));
253 vtkIdTypeArray *d5(vtkIdTypeArray::SafeDownCast(data));
254 if(d2 || d3 || d4 || d5)
255 return ConvertVTKArrayToMCArrayInt(data);
256 std::ostringstream oss;
257 oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
258 throw MZCException(oss.str());
261 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
263 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
264 subMesh->setCoords(coords); subMesh->allocateCells();
265 int nbCells(ca->GetNumberOfCells());
268 //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
269 const vtkIdType *conn(ca->GetData()->GetPointer(0));
270 for(int i=0;i<nbCells;i++)
272 mcIdType sz(ToIdType(*conn++));
273 std::vector<mcIdType> conn2(sz);
274 for(int jj=0;jj<sz;jj++)
275 conn2[jj]=ToIdType(conn[jj]);
276 subMesh->insertNextCell(type,sz,&conn2[0]);
279 return subMesh.retn();
282 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayIdType>& ids)
284 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
285 subMesh->setCoords(coords); subMesh->allocateCells();
286 int nbCells(ca->GetNumberOfCells());
289 //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
290 const vtkIdType *conn(ca->GetData()->GetPointer(0));
291 ids=DataArrayIdType::New() ; ids->alloc(0,1);
292 for(int i=0;i<nbCells;i++)
298 for(int j=0;j<nbTri;j++,conn++)
300 mcIdType conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
301 subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
302 ids->pushBackSilent(i);
307 std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
308 throw MZCException(oss.str());
312 return subMesh.retn();
318 MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
319 MicroField(const std::vector< MicroField >& vs);
320 void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
321 MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
322 std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
324 MCAuto<MEDCouplingUMesh> _m;
325 std::vector<MCAuto<DataArray> > _cellFs;
326 std::vector<MCAuto<DataArray> > _nodeFs;
329 MicroField::MicroField(const std::vector< MicroField >& vs)
331 std::size_t sz(vs.size());
332 std::vector<const MEDCouplingUMesh *> vs2(sz);
333 std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
335 for(std::size_t ii=0;ii<sz;ii++)
337 vs2[ii]=vs[ii].getMesh();
338 arrs2[ii]=vs[ii].getCellFields();
340 nbElts=(int)arrs2[ii].size();
342 if((int)arrs2[ii].size()!=nbElts)
343 throw MZCException("MicroField cstor : internal error !");
345 _cellFs.resize(nbElts);
346 for(int ii=0;ii<nbElts;ii++)
348 std::vector<const DataArray *> arrsTmp(sz);
349 for(std::size_t jj=0;jj<sz;jj++)
351 arrsTmp[jj]=arrs2[jj][ii];
353 _cellFs[ii]=DataArray::Aggregate(arrsTmp);
355 _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
359 void AppendToFields(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, const DataArrayIdType *n2oPtr, typename MEDFileVTKTraits<T>::MCType *dadPtr, MEDFileFields *fs, double timeStep, int tsId)
361 std::string fieldName(dadPtr->getName());
362 MCAuto< typename Traits<T>::FieldType > f(Traits<T>::FieldType::New(tf));
363 f->setTime(timeStep,tsId,0);
365 std::string fieldNameForChuckNorris(MEDCoupling::MEDFileAnyTypeField1TSWithoutSDA::FieldNameToMEDFileConvention(fieldName));
366 f->setName(fieldNameForChuckNorris);
372 MCAuto< typename Traits<T>::ArrayType > dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
376 MCAuto< typename MLFieldTraits<T>::FMTSType > fmts(MLFieldTraits<T>::FMTSType::New());
377 MCAuto< typename MLFieldTraits<T>::F1TSType > f1ts(MLFieldTraits<T>::F1TSType::New());
378 f1ts->setFieldNoProfileSBT(f);
379 fmts->pushBackTimeStep(f1ts);
383 void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayIdType *n2oPtr, double timeStep, int tsId)
385 static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
386 static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
387 if(!da || !mesh || !mfd)
388 throw MZCException("AppendMCFieldFrom : internal error !");
389 MEDFileFields *fs(mfd->getFields());
390 MEDFileMeshes *ms(mfd->getMeshes());
392 throw MZCException("AppendMCFieldFrom : internal error 2 !");
393 MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
396 AppendToFields<double>(tf,mesh,n2oPtr,dad,fs,timeStep,tsId);
399 MCAuto<DataArrayFloat> daf(MEDCoupling::DynamicCast<DataArray,DataArrayFloat>(da));
402 AppendToFields<float>(tf,mesh,n2oPtr,daf,fs,timeStep,tsId);
405 MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
406 MCAuto<DataArrayIdType> daId(MEDCoupling::DynamicCast<DataArray,DataArrayIdType>(da));
407 if(dai.isNotNull() || daId.isNotNull())
409 std::string fieldName(da->getName());
410 if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
413 AppendToFields<mcIdType>(tf,mesh,n2oPtr,daId,fs,timeStep,tsId);
415 AppendToFields<int>(tf,mesh,n2oPtr,dai,fs,timeStep,tsId);
418 else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
420 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
422 throw MZCException("AppendMCFieldFrom : internal error 3 !");
424 throw MZCException("AppendMCFieldFrom : internal error 3 (not mcIdType) !");
426 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daId);
429 MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
430 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
433 else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
435 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
437 throw MZCException("AppendMCFieldFrom : internal error 4 !");
439 throw MZCException("AppendMCFieldFrom : internal error 4 (not mcIdType) !");
441 mm->setFamilyFieldArr(1,daId);
444 MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
445 mm->setFamilyFieldArr(1,dai2);
451 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf, double timeStep, int tsId)
454 throw MZCException("PutAtLevelDealOrder : internal error !");
455 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
456 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
458 throw MZCException("PutAtLevelDealOrder : internal error 2 !");
459 MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
460 mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
461 MCAuto<DataArrayIdType> o2n(mesh->sortCellsInMEDFileFrmt());
462 //const DataArrayIdType *o2nPtr(o2n); // todo: unused
463 MCAuto<DataArrayIdType> n2o;
464 mmu->setMeshAtLevel(meshDimRel,mesh);
465 const DataArrayIdType *n2oPtr(0);
468 n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
470 if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
473 mm->setRenumFieldArr(meshDimRel,n2o);
476 std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
477 for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
479 MCAuto<DataArray> da(*it);
480 AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr,timeStep,tsId);
484 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms, double timeStep, int tsId)
487 throw MZCException("AssignSingleGTMeshes : internal error !");
488 MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
489 MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
491 throw MZCException("AssignSingleGTMeshes : internal error 2 !");
492 int meshDim(-std::numeric_limits<int>::max());
493 std::map<int, std::vector< MicroField > > ms2;
494 for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
496 const MEDCouplingUMesh *elt((*it).getMesh());
499 int myMeshDim(elt->getMeshDimension());
500 meshDim=std::max(meshDim,myMeshDim);
501 ms2[myMeshDim].push_back(*it);
506 for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
508 const std::vector< MicroField >& vs((*it).second);
511 PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0],timeStep,tsId);
515 MicroField merge(vs);
516 PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge,timeStep,tsId);
521 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
524 throw MZCException("BuildCoordsFrom : internal error !");
525 vtkPoints *pts(ds->GetPoints());
527 throw MZCException("BuildCoordsFrom : internal error 2 !");
528 vtkDataArray *data(pts->GetData());
530 throw MZCException("BuildCoordsFrom : internal error 3 !");
531 return ConvertVTKArrayToMCArrayDoubleForced(data);
534 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa, double timeStep, int tsId)
537 throw MZCException("AddNodeFields : internal error !");
538 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
539 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
541 throw MZCException("AddNodeFields : internal error 2 !");
542 MCAuto<MEDCouplingUMesh> mesh;
543 if(!mmu->getNonEmptyLevels().empty())
544 mesh=mmu->getMeshAtLevel(0);
547 mesh=MEDCouplingUMesh::Build0DMeshFromCoords(mmu->getCoords());
548 mesh->setName(mmu->getName());
550 int nba(dsa->GetNumberOfArrays());
551 for(int i=0;i<nba;i++)
553 vtkDataArray *arr(dsa->GetArray(i));
554 const char *name(arr->GetName());
557 MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
559 AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,NULL,timeStep,tsId);
563 std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayIdType *part, vtkDataSetAttributes *dsa)
565 std::vector< MCAuto<DataArray> > ret;
568 int nba(dsa->GetNumberOfArrays());
569 for(int i=0;i<nba;i++)
571 vtkDataArray *arr(dsa->GetArray(i));
574 const char *name(arr->GetName());
575 //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
576 //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
577 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
579 mcarr=mcarr->selectByTupleId(part->begin(),part->end());
580 mcarr->setName(name);
581 ret.push_back(mcarr);
586 std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
588 std::vector< MCAuto<DataArray> > ret;
591 int nba(dsa->GetNumberOfArrays());
592 for(int i=0;i<nba;i++)
594 vtkDataArray *arr(dsa->GetArray(i));
597 const char *name(arr->GetName());
598 //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
599 //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
600 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
601 mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
602 mcarr->setName(name);
603 ret.push_back(mcarr);
608 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
611 throw MZCException("ConvertFromRectilinearGrid : internal error !");
613 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
614 ret->setMeshes(meshes);
615 MCAuto<MEDFileFields> fields(MEDFileFields::New());
616 ret->setFields(fields);
618 MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
619 meshes->pushMesh(cmesh);
620 MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
621 vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
624 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cx));
625 cmeshmc->setCoordsAt(0,arr);
629 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cy));
630 cmeshmc->setCoordsAt(1,arr);
634 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cz));
635 cmeshmc->setCoordsAt(2,arr);
637 std::string meshName(GetMeshNameWithContext(context));
638 cmeshmc->setName(meshName);
639 cmesh->setMesh(cmeshmc);
640 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
641 for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
643 MCAuto<DataArray> da(*it);
644 AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,NULL,timeStep,tsId);
646 std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
647 for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
649 MCAuto<DataArray> da(*it);
650 AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,NULL,timeStep,tsId);
654 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context, double timeStep, int tsId)
657 throw MZCException("ConvertFromPolyData : internal error !");
659 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
660 ret->setMeshes(meshes);
661 MCAuto<MEDFileFields> fields(MEDFileFields::New());
662 ret->setFields(fields);
664 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
665 meshes->pushMesh(umesh);
666 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
667 umesh->setCoords(coords);
668 umesh->setName(GetMeshNameWithContext(context));
671 std::vector< MicroField > ms;
672 vtkCellArray *cd(ds->GetVerts());
675 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
676 if((const MEDCouplingUMesh *)subMesh)
678 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
679 offset+=subMesh->getNumberOfCells();
680 ms.push_back(MicroField(subMesh,cellFs));
683 vtkCellArray *cc(ds->GetLines());
686 MCAuto<MEDCouplingUMesh> subMesh;
689 subMesh=BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2);
691 catch(INTERP_KERNEL::Exception& e)
693 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;
694 throw INTERP_KERNEL::Exception(oss.str());
696 if((const MEDCouplingUMesh *)subMesh)
698 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
699 offset+=subMesh->getNumberOfCells();
700 ms.push_back(MicroField(subMesh,cellFs));
703 vtkCellArray *cb(ds->GetPolys());
706 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
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 *ca(ds->GetStrips());
717 MCAuto<DataArrayIdType> ids;
718 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
719 if((const MEDCouplingUMesh *)subMesh)
721 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
722 offset+=subMesh->getNumberOfCells();
723 ms.push_back(MicroField(subMesh,cellFs));
726 AssignSingleGTMeshes(ret,ms,timeStep,tsId);
727 AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
730 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
733 throw MZCException("ConvertFromUnstructuredGrid : internal error !");
735 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
736 ret->setMeshes(meshes);
737 MCAuto<MEDFileFields> fields(MEDFileFields::New());
738 ret->setFields(fields);
740 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
741 meshes->pushMesh(umesh);
742 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
743 umesh->setCoords(coords);
744 umesh->setName(GetMeshNameWithContext(context));
745 vtkIdType nbCells(ds->GetNumberOfCells());
746 vtkCellArray *ca(ds->GetCells());
749 //vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries()); // todo: unused
750 //vtkIdType *caPtr(ca->GetData()->GetPointer(0)); // todo: unused
751 vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
753 throw MZCException("ConvertFromUnstructuredGrid : internal error");
754 vtkIdTypeArray *cla(ds->GetCellLocationsArray());
755 //const vtkIdType *claPtr(cla->GetPointer(0)); // todo: unused
757 throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
758 const unsigned char *ctPtr(ct->GetPointer(0));
759 std::map<int,int> m(ComputeMapOfType());
760 MCAuto<DataArrayInt> lev(DataArrayInt::New()) ; lev->alloc(nbCells,1);
761 int *levPtr(lev->getPointer());
762 for(vtkIdType i=0;i<nbCells;i++)
764 std::map<int,int>::iterator it(m.find(ctPtr[i]));
767 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
768 levPtr[i]=cm.getDimension();
772 if(ctPtr[i]==VTK_POLY_VERTEX)
774 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(INTERP_KERNEL::NORM_POINT1));
775 levPtr[i]=cm.getDimension();
779 std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
780 throw MZCException(oss.str());
784 MCAuto<DataArrayInt> levs(lev->getDifferentValues());
785 std::vector< MicroField > ms;
786 //vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations()); // todo: unused
787 for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
789 MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
790 m0->setCoords(coords); m0->allocateCells();
791 MCAuto<DataArrayIdType> cellIdsCurLev(lev->findIdsEqual(*curLev));
792 for(const mcIdType *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
794 int vtkType(ds->GetCellType(*cellId));
795 std::map<int,int>::iterator it(m.find(vtkType));
796 INTERP_KERNEL::NormalizedCellType ct=it!=m.end()?(INTERP_KERNEL::NormalizedCellType)((*it).second):INTERP_KERNEL::NORM_POINT1;
797 if(ct!=INTERP_KERNEL::NORM_POLYHED && vtkType!=VTK_POLY_VERTEX)
800 const vtkIdType *pts(nullptr);
801 ds->GetCellPoints(*cellId, sz, pts);
802 std::vector<mcIdType> conn2(pts,pts+sz);
803 m0->insertNextCell(ct,sz,conn2.data());
805 else if(ct==INTERP_KERNEL::NORM_POLYHED)
807 // # de faces du polyèdre
808 vtkIdType nbOfFaces(0);
809 // connectivé des faces (numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...)
810 const vtkIdType *facPtr(nullptr);
811 ds->GetFaceStream(*cellId, nbOfFaces, facPtr);
812 std::vector<mcIdType> conn;
813 for(vtkIdType k=0;k<nbOfFaces;k++)
815 vtkIdType nbOfNodesInFace(*facPtr++);
816 std::copy(facPtr,facPtr+nbOfNodesInFace,std::back_inserter(conn));
819 facPtr+=nbOfNodesInFace;
821 m0->insertNextCell(ct,ToIdType(conn.size()),&conn[0]);
826 const vtkIdType *pts(nullptr);
827 ds->GetCellPoints(*cellId, sz, pts);
829 throw MZCException("ConvertFromUnstructuredGrid : non single poly vertex not managed by MED !");
830 m0->insertNextCell(ct,1,(const mcIdType*)pts);
833 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
834 ms.push_back(MicroField(m0,cellFs));
836 AssignSingleGTMeshes(ret,ms,timeStep,tsId);
837 AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
842 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
845 throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
846 vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
847 vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
848 vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
851 ConvertFromPolyData(mfd,ds2,context,timeStep,tsId);
855 ConvertFromUnstructuredGrid(mfd,ds3,context,timeStep,tsId);
859 ConvertFromRectilinearGrid(mfd,ds4,context,timeStep,tsId);
862 throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
865 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
868 throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
869 int nbBlocks(ds->GetNumberOfBlocks());
870 if(nbBlocks==1 && context.empty())
872 vtkDataObject *uniqueElt(ds->GetBlock(0));
874 throw MZCException("Unique elt in multiblock is NULL !");
875 vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
878 WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context,timeStep,tsId);
882 for(int i=0;i<nbBlocks;i++)
884 vtkDataObject *elt(ds->GetBlock(i));
885 std::vector<int> context2;
886 context2.push_back(i);
889 std::ostringstream oss; oss << "In context ";
890 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
891 oss << " at pos #" << i << " elt is NULL !";
892 throw MZCException(oss.str());
894 vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
897 WriteMEDFileFromVTKDataSet(mfd,elt1,context,timeStep,tsId);
900 vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
903 WriteMEDFileFromVTKMultiBlock(mfd,elt2,context,timeStep,tsId);
906 std::ostringstream oss; oss << "In context ";
907 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
908 oss << " at pos #" << i << " elt not recognized data type !";
909 throw MZCException(oss.str());
913 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input, double timeStep, int tsId)
916 throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
917 std::vector<int> context;
918 vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
921 WriteMEDFileFromVTKDataSet(mfd,input1,context,timeStep,tsId);
924 vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
927 WriteMEDFileFromVTKMultiBlock(mfd,input2,context,timeStep,tsId);
930 throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
933 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
939 MEDFileMeshes *meshes(mfd->getMeshes());
942 if(meshes->getNumberOfMeshes()!=1)
944 MEDFileMesh *mm(meshes->getMeshAtPos(0));
947 mm->setName(meshName);
948 for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
949 mm->setFamilyId((*it).getName(),(*it).getID());
950 for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
951 mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
952 MEDFileFields *fields(mfd->getFields());
955 for(int i=0;i<fields->getNumberOfFields();i++)
957 MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
960 fmts->setMeshName(meshName);