1 // Copyright (C) 2016 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay (EDF R&D)
21 #include "vtkMEDWriter.h"
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkIntArray.h"
25 #include "vtkCellData.h"
26 #include "vtkPointData.h"
27 #include "vtkFloatArray.h"
28 #include "vtkCellArray.h"
30 #include "vtkStreamingDemandDrivenPipeline.h"
31 #include "vtkInformationDataObjectMetaDataKey.h"
32 #include "vtkUnstructuredGrid.h"
33 #include "vtkMultiBlockDataSet.h"
34 #include "vtkRectilinearGrid.h"
35 #include "vtkInformationStringKey.h"
36 #include "vtkAlgorithmOutput.h"
37 #include "vtkObjectFactory.h"
38 #include "vtkMutableDirectedGraph.h"
39 #include "vtkMultiBlockDataSet.h"
40 #include "vtkPolyData.h"
41 #include "vtkDataSet.h"
42 #include "vtkInformationVector.h"
43 #include "vtkInformation.h"
44 #include "vtkDataArraySelection.h"
45 #include "vtkTimeStamp.h"
46 #include "vtkInEdgeIterator.h"
47 #include "vtkInformationDataObjectKey.h"
48 #include "vtkExecutive.h"
49 #include "vtkVariantArray.h"
50 #include "vtkStringArray.h"
51 #include "vtkDoubleArray.h"
52 #include "vtkCharArray.h"
53 #include "vtkUnsignedCharArray.h"
54 #include "vtkDataSetAttributes.h"
55 #include "vtkDemandDrivenPipeline.h"
56 #include "vtkDataObjectTreeIterator.h"
57 #include "vtkWarpScalar.h"
59 #include "MEDFileMesh.hxx"
60 #include "MEDFileField.hxx"
61 #include "MEDFileData.hxx"
62 #include "MEDCouplingMemArray.hxx"
63 #include "MEDCouplingFieldInt.hxx"
64 #include "MEDCouplingFieldDouble.hxx"
65 #include "MEDCouplingRefCountObject.hxx"
71 using MEDCoupling::MEDFileData;
72 using MEDCoupling::MEDFileMesh;
73 using MEDCoupling::MEDFileCMesh;
74 using MEDCoupling::MEDFileUMesh;
75 using MEDCoupling::MEDFileFields;
76 using MEDCoupling::MEDFileMeshes;
78 using MEDCoupling::MEDFileIntField1TS;
79 using MEDCoupling::MEDFileField1TS;
80 using MEDCoupling::MEDFileIntFieldMultiTS;
81 using MEDCoupling::MEDFileFieldMultiTS;
82 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
83 using MEDCoupling::DataArray;
84 using MEDCoupling::DataArrayInt;
85 using MEDCoupling::DataArrayDouble;
86 using MEDCoupling::MEDCouplingMesh;
87 using MEDCoupling::MEDCouplingUMesh;
88 using MEDCoupling::MEDCouplingCMesh;
89 using MEDCoupling::MEDCouplingFieldDouble;
90 using MEDCoupling::MEDCouplingFieldInt;
91 using MEDCoupling::MCAuto;
93 vtkStandardNewMacro(vtkMEDWriter);
97 class MZCException : public std::exception
100 MZCException(const std::string& s):_reason(s) { }
101 virtual const char *what() const throw() { return _reason.c_str(); }
102 virtual ~MZCException() throw() { }
109 std::map<int,int> ComputeMapOfType()
111 std::map<int,int> ret;
112 int nbOfTypesInMC(sizeof(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER)/sizeof(int));
113 for(int i=0;i<nbOfTypesInMC;i++)
115 int vtkId(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER[i]);
122 std::string GetMeshNameWithContext(const std::vector<int>& context)
124 static const char DFT_MESH_NAME[]="Mesh";
126 return DFT_MESH_NAME;
127 std::ostringstream oss; oss << DFT_MESH_NAME;
128 for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
133 DataArrayInt *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
136 throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
137 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
138 std::size_t nbElts(nbTuples*nbComp);
139 MCAuto<DataArrayInt> ret(DataArrayInt::New());
140 ret->alloc(nbTuples,nbComp);
141 for(int i=0;i<nbComp;i++)
143 const char *comp(data->GetComponentName(i));
145 ret->setInfoOnComponent(i,comp);
147 int *ptOut(ret->getPointer());
148 vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
151 const int *pt(d0->GetPointer(0));
152 std::copy(pt,pt+nbElts,ptOut);
155 std::ostringstream oss;
156 oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
157 throw MZCException(oss.str());
160 DataArrayDouble *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
163 throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
164 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
165 std::size_t nbElts(nbTuples*nbComp);
166 MCAuto<DataArrayDouble> ret(DataArrayDouble::New());
167 ret->alloc(nbTuples,nbComp);
168 for(int i=0;i<nbComp;i++)
170 const char *comp(data->GetComponentName(i));
172 ret->setInfoOnComponent(i,comp);
174 double *ptOut(ret->getPointer());
175 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
178 const float *pt(d0->GetPointer(0));
179 for(std::size_t i=0;i<nbElts;i++)
183 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
186 const double *pt(d1->GetPointer(0));
187 std::copy(pt,pt+nbElts,ptOut);
190 std::ostringstream oss;
191 oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << typeid(*data).name() << "\" type !";
192 throw MZCException(oss.str());
195 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
198 throw MZCException("ConvertVTKArrayToMCArray : internal error !");
199 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
200 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
202 return ConvertVTKArrayToMCArrayDouble(data);
203 vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
205 return ConvertVTKArrayToMCArrayInt(data);
206 std::ostringstream oss;
207 oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
208 throw MZCException(oss.str());
211 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
213 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
214 subMesh->setCoords(coords); subMesh->allocateCells();
215 int nbCells(ca->GetNumberOfCells());
218 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
219 const vtkIdType *conn(ca->GetPointer());
220 for(int i=0;i<nbCells;i++)
223 std::vector<int> conn2(sz);
224 for(int jj=0;jj<sz;jj++)
226 subMesh->insertNextCell(type,sz,&conn2[0]);
229 return subMesh.retn();
232 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayInt>& ids)
234 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
235 subMesh->setCoords(coords); subMesh->allocateCells();
236 int nbCells(ca->GetNumberOfCells());
239 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
240 const vtkIdType *conn(ca->GetPointer());
241 ids=DataArrayInt::New() ; ids->alloc(0,1);
242 for(int i=0;i<nbCells;i++)
248 for(int j=0;j<nbTri;j++,conn++)
250 int conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
251 subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
252 ids->pushBackSilent(i);
257 std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
258 throw MZCException(oss.str());
262 return subMesh.retn();
268 MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
269 MicroField(const std::vector< MicroField >& vs);
270 void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
271 MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
272 std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
274 MCAuto<MEDCouplingUMesh> _m;
275 std::vector<MCAuto<DataArray> > _cellFs;
276 std::vector<MCAuto<DataArray> > _nodeFs;
279 MicroField::MicroField(const std::vector< MicroField >& vs)
281 std::size_t sz(vs.size());
282 std::vector<const MEDCouplingUMesh *> vs2(sz);
283 std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
285 for(std::size_t ii=0;ii<sz;ii++)
287 vs2[ii]=vs[ii].getMesh();
288 arrs2[ii]=vs[ii].getCellFields();
290 nbElts=arrs2[ii].size();
292 if(arrs2[ii].size()!=nbElts)
293 throw MZCException("MicroField cstor : internal error !");
295 _cellFs.resize(nbElts);
296 for(int ii=0;ii<nbElts;ii++)
298 std::vector<const DataArray *> arrsTmp(sz);
299 for(std::size_t jj=0;jj<sz;jj++)
301 arrsTmp[jj]=arrs2[jj][ii];
303 _cellFs[ii]=DataArray::Aggregate(arrsTmp);
305 _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
308 void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayInt *n2oPtr)
310 static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
311 static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
312 if(!da || !mesh || !mfd)
313 throw MZCException("AppendMCFieldFrom : internal error !");
314 MEDFileFields *fs(mfd->getFields());
315 MEDFileMeshes *ms(mfd->getMeshes());
317 throw MZCException("AppendMCFieldFrom : internal error 2 !");
318 MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
319 DataArrayDouble *dadPtr(dad);
320 std::string fieldName;
323 fieldName=dadPtr->getName();
324 MCAuto<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
325 f->setName(fieldName);
330 MCAuto<DataArrayDouble> dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
334 MCAuto<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
335 MCAuto<MEDFileField1TS> f1ts(MEDFileField1TS::New());
336 f1ts->setFieldNoProfileSBT(f);
337 fmts->pushBackTimeStep(f1ts);
341 MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
342 DataArrayInt *daiPtr(dai);
345 fieldName=daiPtr->getName();
346 if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
348 MCAuto<MEDCouplingFieldInt> f(MEDCouplingFieldInt::New(tf));
349 f->setName(fieldName);
351 MCAuto<MEDFileIntFieldMultiTS> fmts(MEDFileIntFieldMultiTS::New());
352 MCAuto<MEDFileIntField1TS> f1ts(MEDFileIntField1TS::New());
356 f1ts->setFieldNoProfileSBT(f);
360 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
362 f1ts->setFieldNoProfileSBT(f);
364 fmts->pushBackTimeStep(f1ts);
368 else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
370 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
372 throw MZCException("AppendMCFieldFrom : internal error 3 !");
374 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daiPtr);
377 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
378 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
381 else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
383 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
385 throw MZCException("AppendMCFieldFrom : internal error 4 !");
387 mm->setFamilyFieldArr(1,daiPtr);
390 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
391 mm->setFamilyFieldArr(1,dai2);
397 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf)
400 throw MZCException("PutAtLevelDealOrder : internal error !");
401 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
402 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
404 throw MZCException("PutAtLevelDealOrder : internal error 2 !");
405 MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
406 mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
407 MCAuto<DataArrayInt> o2n(mesh->sortCellsInMEDFileFrmt());
408 const DataArrayInt *o2nPtr(o2n);
409 MCAuto<DataArrayInt> n2o;
410 mmu->setMeshAtLevel(meshDimRel,mesh);
411 const DataArrayInt *n2oPtr(0);
414 n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
416 if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
419 mm->setRenumFieldArr(meshDimRel,n2o);
422 std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
423 for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
425 MCAuto<DataArray> da(*it);
426 AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr);
430 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms)
433 throw MZCException("AssignSingleGTMeshes : internal error !");
434 MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
435 MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
437 throw MZCException("AssignSingleGTMeshes : internal error 2 !");
438 int meshDim(-std::numeric_limits<int>::max());
439 std::map<int, std::vector< MicroField > > ms2;
440 for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
442 const MEDCouplingUMesh *elt((*it).getMesh());
445 int myMeshDim(elt->getMeshDimension());
446 meshDim=std::max(meshDim,myMeshDim);
447 ms2[myMeshDim].push_back(*it);
452 for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
454 const std::vector< MicroField >& vs((*it).second);
457 PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0]);
461 MicroField merge(vs);
462 PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge);
467 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
470 throw MZCException("BuildCoordsFrom : internal error !");
471 vtkPoints *pts(ds->GetPoints());
473 throw MZCException("BuildCoordsFrom : internal error 2 !");
474 vtkDataArray *data(pts->GetData());
476 throw MZCException("BuildCoordsFrom : internal error 3 !");
477 MCAuto<DataArrayDouble> coords(ConvertVTKArrayToMCArrayDouble(data));
478 return coords.retn();
481 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa)
484 throw MZCException("AddNodeFields : internal error !");
485 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
486 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
488 throw MZCException("AddNodeFields : internal error 2 !");
489 MCAuto<MEDCouplingUMesh> mesh(mmu->getMeshAtLevel(0));
490 int nba(dsa->GetNumberOfArrays());
491 for(int i=0;i<nba;i++)
493 vtkDataArray *arr(dsa->GetArray(i));
494 const char *name(arr->GetName());
497 MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
499 AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,0);
503 std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayInt *part, vtkDataSetAttributes *dsa)
505 std::vector< MCAuto<DataArray> > ret;
508 int nba(dsa->GetNumberOfArrays());
509 for(int i=0;i<nba;i++)
511 vtkDataArray *arr(dsa->GetArray(i));
514 const char *name(arr->GetName());
515 int nbCompo(arr->GetNumberOfComponents());
516 vtkIdType nbTuples(arr->GetNumberOfTuples());
517 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
519 mcarr=mcarr->selectByTupleId(part->begin(),part->end());
520 mcarr->setName(name);
521 ret.push_back(mcarr);
526 std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
528 std::vector< MCAuto<DataArray> > ret;
531 int nba(dsa->GetNumberOfArrays());
532 for(int i=0;i<nba;i++)
534 vtkDataArray *arr(dsa->GetArray(i));
537 const char *name(arr->GetName());
538 int nbCompo(arr->GetNumberOfComponents());
539 vtkIdType nbTuples(arr->GetNumberOfTuples());
540 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
541 mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
542 mcarr->setName(name);
543 ret.push_back(mcarr);
548 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context)
551 throw MZCException("ConvertFromRectilinearGrid : internal error !");
553 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
554 ret->setMeshes(meshes);
555 MCAuto<MEDFileFields> fields(MEDFileFields::New());
556 ret->setFields(fields);
558 MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
559 meshes->pushMesh(cmesh);
560 MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
561 vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
564 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cx));
565 cmeshmc->setCoordsAt(0,arr);
569 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cy));
570 cmeshmc->setCoordsAt(1,arr);
574 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cz));
575 cmeshmc->setCoordsAt(2,arr);
577 std::string meshName(GetMeshNameWithContext(context));
578 cmeshmc->setName(meshName);
579 cmesh->setMesh(cmeshmc);
580 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
581 for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
583 MCAuto<DataArray> da(*it);
584 AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,0);
586 std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
587 for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
589 MCAuto<DataArray> da(*it);
590 AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,0);
594 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context)
597 throw MZCException("ConvertFromPolyData : internal error !");
599 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
600 ret->setMeshes(meshes);
601 MCAuto<MEDFileFields> fields(MEDFileFields::New());
602 ret->setFields(fields);
604 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
605 meshes->pushMesh(umesh);
606 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
607 umesh->setCoords(coords);
608 umesh->setName(GetMeshNameWithContext(context));
611 std::vector< MicroField > ms;
612 vtkCellArray *cd(ds->GetVerts());
615 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
616 if((const MEDCouplingUMesh *)subMesh)
618 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
619 offset+=subMesh->getNumberOfCells();
620 ms.push_back(MicroField(subMesh,cellFs));
623 vtkCellArray *cc(ds->GetLines());
626 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2));
627 if((const MEDCouplingUMesh *)subMesh)
629 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
630 offset+=subMesh->getNumberOfCells();
631 ms.push_back(MicroField(subMesh,cellFs));
634 vtkCellArray *cb(ds->GetPolys());
637 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
638 if((const MEDCouplingUMesh *)subMesh)
640 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
641 offset+=subMesh->getNumberOfCells();
642 ms.push_back(MicroField(subMesh,cellFs));
645 vtkCellArray *ca(ds->GetStrips());
648 MCAuto<DataArrayInt> ids;
649 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
650 if((const MEDCouplingUMesh *)subMesh)
652 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
653 offset+=subMesh->getNumberOfCells();
654 ms.push_back(MicroField(subMesh,cellFs));
657 AssignSingleGTMeshes(ret,ms);
658 AddNodeFields(ret,ds->GetPointData());
661 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context)
664 throw MZCException("ConvertFromUnstructuredGrid : internal error !");
666 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
667 ret->setMeshes(meshes);
668 MCAuto<MEDFileFields> fields(MEDFileFields::New());
669 ret->setFields(fields);
671 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
672 meshes->pushMesh(umesh);
673 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
674 umesh->setCoords(coords);
675 umesh->setName(GetMeshNameWithContext(context));
676 vtkIdType nbCells(ds->GetNumberOfCells());
677 vtkCellArray *ca(ds->GetCells());
680 vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries());
681 vtkIdType *caPtr(ca->GetPointer());
682 vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
684 throw MZCException("ConvertFromUnstructuredGrid : internal error");
685 vtkIdTypeArray *cla(ds->GetCellLocationsArray());
686 const vtkIdType *claPtr(cla->GetPointer(0));
688 throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
689 const unsigned char *ctPtr(ct->GetPointer(0));
690 std::map<int,int> m(ComputeMapOfType());
691 MCAuto<DataArrayInt> lev(DataArrayInt::New()) ; lev->alloc(nbCells,1);
692 int *levPtr(lev->getPointer());
693 for(vtkIdType i=0;i<nbCells;i++)
695 std::map<int,int>::iterator it(m.find(ctPtr[i]));
698 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
699 levPtr[i]=cm.getDimension();
703 std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
704 throw MZCException(oss.str());
708 int meshDim(lev->getMaxValue(dummy));
709 MCAuto<DataArrayInt> levs(lev->getDifferentValues());
710 std::vector< MicroField > ms;
711 vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations());
712 for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
714 MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
715 m0->setCoords(coords); m0->allocateCells();
716 MCAuto<DataArrayInt> cellIdsCurLev(lev->findIdsEqual(*curLev));
717 for(const int *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
719 std::map<int,int>::iterator it(m.find(ctPtr[*cellId]));
720 vtkIdType offset(claPtr[*cellId]);
721 vtkIdType sz(caPtr[offset]);
722 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)(*it).second);
723 if(ct!=INTERP_KERNEL::NORM_POLYHED)
725 std::vector<int> conn2(sz);
726 for(int kk=0;kk<sz;kk++)
727 conn2[kk]=caPtr[offset+1+kk];
728 m0->insertNextCell(ct,sz,&conn2[0]);
732 if(!faces || !faceLoc)
733 throw MZCException("ConvertFromUnstructuredGrid : faces are expected when there are polyhedra !");
734 const vtkIdType *facPtr(faces->GetPointer(0)),*facLocPtr(faceLoc->GetPointer(0));
735 std::vector<int> conn;
736 int off0(facLocPtr[*cellId]);
737 int nbOfFaces(facPtr[off0++]);
738 for(int k=0;k<nbOfFaces;k++)
740 int nbOfNodesInFace(facPtr[off0++]);
741 std::copy(facPtr+off0,facPtr+off0+nbOfNodesInFace,std::back_inserter(conn));
742 off0+=nbOfNodesInFace;
746 m0->insertNextCell(ct,conn.size(),&conn[0]);
749 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
750 ms.push_back(MicroField(m0,cellFs));
752 AssignSingleGTMeshes(ret,ms);
753 AddNodeFields(ret,ds->GetPointData());
758 vtkMEDWriter::vtkMEDWriter():WriteAllTimeSteps(0),FileName(0),IsTouched(false)
762 vtkMEDWriter::~vtkMEDWriter()
766 vtkInformationDataObjectMetaDataKey *GetMEDReaderMetaDataIfAny()
768 static const char ZE_KEY[]="vtkMEDReader::META_DATA";
769 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
770 if(!gd->hasKey(ZE_KEY))
772 std::string ptSt(gd->value(ZE_KEY));
774 std::istringstream iss(ptSt); iss >> pt;
775 return reinterpret_cast<vtkInformationDataObjectMetaDataKey *>(pt);
778 bool IsInformationOK(vtkInformation *info)
780 vtkInformationDataObjectMetaDataKey *key(GetMEDReaderMetaDataIfAny());
783 // Check the information contain meta data key
787 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
791 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
792 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
795 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
797 vtkStdString &st(verticesNames2->GetValue(i));
798 if(st=="MeshesFamsGrps")
807 Grp(const std::string& name):_name(name) { }
808 void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
809 std::string getName() const { return _name; }
810 std::vector<std::string> getFamilies() const { return _fams; }
813 std::vector<std::string> _fams;
819 Fam(const std::string& name)
821 static const char ZE_SEP[]="@@][@@";
822 std::size_t pos(name.find(ZE_SEP));
823 std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
824 std::istringstream iss(name1);
828 std::string getName() const { return _name; }
829 int getID() const { return _id; }
835 void LoadFamGrpMapInfo(vtkMutableDirectedGraph *sil, std::string& meshName, std::vector<Grp>& groups, std::vector<Fam>& fams)
838 throw MZCException("LoadFamGrpMapInfo : internal error !");
840 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
841 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
844 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
846 vtkStdString &st(verticesNames2->GetValue(i));
847 if(st=="MeshesFamsGrps")
854 throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
855 vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
856 sil->GetAdjacentVertices(id0,it0);
858 while(it0->HasNext())
860 vtkIdType id1(it0->Next());
861 std::string mName((const char *)verticesNames2->GetValue(id1));
863 vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
864 sil->GetAdjacentVertices(id1,it1);
865 vtkIdType idZeGrps(it1->Next());//zeGroups
866 vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
867 sil->GetAdjacentVertices(idZeGrps,itGrps);
868 while(itGrps->HasNext())
870 vtkIdType idg(itGrps->Next());
871 Grp grp((const char *)verticesNames2->GetValue(idg));
872 vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
873 sil->GetAdjacentVertices(idg,itGrps2);
874 std::vector<std::string> famsOnGroup;
875 while(itGrps2->HasNext())
877 vtkIdType idgf(itGrps2->Next());
878 famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
880 grp.setFamilies(famsOnGroup);
882 groups.push_back(grp);
885 vtkIdType idZeFams(it1->Next());//zeFams
887 vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
888 sil->GetAdjacentVertices(idZeFams,itFams);
889 while(itFams->HasNext())
891 vtkIdType idf(itFams->Next());
892 Fam fam((const char *)verticesNames2->GetValue(idf));
900 int vtkMEDWriter::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
902 //std::cerr << "########################################## vtkMEDWriter::RequestInformation ########################################## " << (const char *) this->FileName << std::endl;
906 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context)
909 throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
910 vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
911 vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
912 vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
915 ConvertFromPolyData(mfd,ds2,context);
919 ConvertFromUnstructuredGrid(mfd,ds3,context);
923 ConvertFromRectilinearGrid(mfd,ds4,context);
926 throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
929 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context)
932 throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
933 int nbBlocks(ds->GetNumberOfBlocks());
934 if(nbBlocks==1 && context.empty())
936 vtkDataObject *uniqueElt(ds->GetBlock(0));
938 throw MZCException("Unique elt in multiblock is NULL !");
939 vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
942 WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context);
946 for(int i=0;i<nbBlocks;i++)
948 vtkDataObject *elt(ds->GetBlock(i));
949 std::vector<int> context2;
950 context2.push_back(i);
953 std::ostringstream oss; oss << "In context ";
954 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
955 oss << " at pos #" << i << " elt is NULL !";
956 throw MZCException(oss.str());
958 vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
961 WriteMEDFileFromVTKDataSet(mfd,elt1,context);
964 vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
967 WriteMEDFileFromVTKMultiBlock(mfd,elt2,context);
970 std::ostringstream oss; oss << "In context ";
971 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
972 oss << " at pos #" << i << " elt not recognized data type !";
973 throw MZCException(oss.str());
977 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input)
980 throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
981 std::vector<int> context;
982 vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
985 WriteMEDFileFromVTKDataSet(mfd,input1,context);
988 vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
991 WriteMEDFileFromVTKMultiBlock(mfd,input2,context);
994 throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
997 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
1001 if(meshName.empty())
1003 MEDFileMeshes *meshes(mfd->getMeshes());
1006 if(meshes->getNumberOfMeshes()!=1)
1008 MEDFileMesh *mm(meshes->getMeshAtPos(0));
1011 mm->setName(meshName);
1012 for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
1013 mm->setFamilyId((*it).getName(),(*it).getID());
1014 for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
1015 mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
1016 MEDFileFields *fields(mfd->getFields());
1019 for(int i=0;i<fields->getNumberOfFields();i++)
1021 MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
1024 fmts->setMeshName(meshName);
1028 int vtkMEDWriter::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
1030 //std::cerr << "########################################## vtkMEDWriter::RequestData ########################################## " << (const char *) this->FileName << std::endl;
1033 vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
1034 std::string meshName;
1035 std::vector<Grp> groups; std::vector<Fam> fams;
1036 if(IsInformationOK(inputInfo))
1038 vtkMutableDirectedGraph *famGrpGraph(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
1039 LoadFamGrpMapInfo(famGrpGraph,meshName,groups,fams);
1041 vtkInformation *outInfo(outputVector->GetInformationObject(0));
1042 vtkDataObject *input(vtkDataObject::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
1044 throw MZCException("Not recognized data object in input of the MEDWriter ! Maybe not implemented yet !");
1045 MCAuto<MEDFileData> mfd(MEDFileData::New());
1046 WriteMEDFileFromVTKGDS(mfd,input);
1047 PutFamGrpInfoIfAny(mfd,meshName,groups,fams);
1048 mfd->write(this->FileName,this->IsTouched?0:2); this->IsTouched=true;
1049 outInfo->Set(vtkDataObject::DATA_OBJECT(),input);
1051 catch(MZCException& e)
1053 std::ostringstream oss;
1054 oss << "Exception has been thrown in vtkMEDWriter::RequestData : During writing of \"" << (const char *) this->FileName << "\", the following exception has been thrown : "<< e.what() << std::endl;
1055 if(this->HasObserver("ErrorEvent") )
1056 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
1058 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
1059 vtkObject::BreakOnError();
1065 void vtkMEDWriter::PrintSelf(ostream& os, vtkIndent indent)
1067 this->Superclass::PrintSelf(os, indent);
1070 int vtkMEDWriter::Write()