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 "MEDCouplingFieldDouble.hxx"
64 #include "MEDCouplingRefCountObject.hxx"
70 using MEDCoupling::MEDFileData;
71 using MEDCoupling::MEDFileMesh;
72 using MEDCoupling::MEDFileCMesh;
73 using MEDCoupling::MEDFileUMesh;
74 using MEDCoupling::MEDFileFields;
75 using MEDCoupling::MEDFileMeshes;
77 using MEDCoupling::MEDFileIntField1TS;
78 using MEDCoupling::MEDFileField1TS;
79 using MEDCoupling::MEDFileIntFieldMultiTS;
80 using MEDCoupling::MEDFileFieldMultiTS;
81 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
82 using MEDCoupling::DataArray;
83 using MEDCoupling::DataArrayInt;
84 using MEDCoupling::DataArrayDouble;
85 using MEDCoupling::MEDCouplingMesh;
86 using MEDCoupling::MEDCouplingUMesh;
87 using MEDCoupling::MEDCouplingCMesh;
88 using MEDCoupling::MEDCouplingFieldDouble;
89 using MEDCoupling::MCAuto;
91 vtkStandardNewMacro(vtkMEDWriter);
95 class MZCException : public std::exception
98 MZCException(const std::string& s):_reason(s) { }
99 virtual const char *what() const throw() { return _reason.c_str(); }
100 virtual ~MZCException() throw() { }
107 std::map<int,int> ComputeMapOfType()
109 std::map<int,int> ret;
110 int nbOfTypesInMC(sizeof(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER)/sizeof(int));
111 for(int i=0;i<nbOfTypesInMC;i++)
113 int vtkId(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER[i]);
120 std::string GetMeshNameWithContext(const std::vector<int>& context)
122 static const char DFT_MESH_NAME[]="Mesh";
124 return DFT_MESH_NAME;
125 std::ostringstream oss; oss << DFT_MESH_NAME;
126 for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
131 DataArrayInt *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
134 throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
135 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
136 std::size_t nbElts(nbTuples*nbComp);
137 MCAuto<DataArrayInt> ret(DataArrayInt::New());
138 ret->alloc(nbTuples,nbComp);
139 for(int i=0;i<nbComp;i++)
141 const char *comp(data->GetComponentName(i));
143 ret->setInfoOnComponent(i,comp);
145 int *ptOut(ret->getPointer());
146 vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
149 const int *pt(d0->GetPointer(0));
150 std::copy(pt,pt+nbElts,ptOut);
153 std::ostringstream oss;
154 oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
155 throw MZCException(oss.str());
158 DataArrayDouble *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
161 throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
162 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
163 std::size_t nbElts(nbTuples*nbComp);
164 MCAuto<DataArrayDouble> ret(DataArrayDouble::New());
165 ret->alloc(nbTuples,nbComp);
166 for(int i=0;i<nbComp;i++)
168 const char *comp(data->GetComponentName(i));
170 ret->setInfoOnComponent(i,comp);
172 double *ptOut(ret->getPointer());
173 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
176 const float *pt(d0->GetPointer(0));
177 for(std::size_t i=0;i<nbElts;i++)
181 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
184 const double *pt(d1->GetPointer(0));
185 std::copy(pt,pt+nbElts,ptOut);
188 std::ostringstream oss;
189 oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << typeid(*data).name() << "\" type !";
190 throw MZCException(oss.str());
193 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
196 throw MZCException("ConvertVTKArrayToMCArray : internal error !");
197 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
198 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
200 return ConvertVTKArrayToMCArrayDouble(data);
201 vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
203 return ConvertVTKArrayToMCArrayInt(data);
204 std::ostringstream oss;
205 oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
206 throw MZCException(oss.str());
209 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
211 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
212 subMesh->setCoords(coords); subMesh->allocateCells();
213 int nbCells(ca->GetNumberOfCells());
216 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
217 const vtkIdType *conn(ca->GetPointer());
218 for(int i=0;i<nbCells;i++)
221 std::vector<int> conn2(sz);
222 for(int jj=0;jj<sz;jj++)
224 subMesh->insertNextCell(type,sz,&conn2[0]);
227 return subMesh.retn();
230 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayInt>& ids)
232 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
233 subMesh->setCoords(coords); subMesh->allocateCells();
234 int nbCells(ca->GetNumberOfCells());
237 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
238 const vtkIdType *conn(ca->GetPointer());
239 ids=DataArrayInt::New() ; ids->alloc(0,1);
240 for(int i=0;i<nbCells;i++)
246 for(int j=0;j<nbTri;j++,conn++)
248 int conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
249 subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
250 ids->pushBackSilent(i);
255 std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
256 throw MZCException(oss.str());
260 return subMesh.retn();
266 MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
267 MicroField(const std::vector< MicroField >& vs);
268 void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
269 MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
270 std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
272 MCAuto<MEDCouplingUMesh> _m;
273 std::vector<MCAuto<DataArray> > _cellFs;
274 std::vector<MCAuto<DataArray> > _nodeFs;
277 MicroField::MicroField(const std::vector< MicroField >& vs)
279 std::size_t sz(vs.size());
280 std::vector<const MEDCouplingUMesh *> vs2(sz);
281 std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
283 for(std::size_t ii=0;ii<sz;ii++)
285 vs2[ii]=vs[ii].getMesh();
286 arrs2[ii]=vs[ii].getCellFields();
288 nbElts=arrs2[ii].size();
290 if(arrs2[ii].size()!=nbElts)
291 throw MZCException("MicroField cstor : internal error !");
293 _cellFs.resize(nbElts);
294 for(int ii=0;ii<nbElts;ii++)
296 std::vector<const DataArray *> arrsTmp(sz);
297 for(std::size_t jj=0;jj<sz;jj++)
299 arrsTmp[jj]=arrs2[jj][ii];
301 _cellFs[ii]=DataArray::Aggregate(arrsTmp);
303 _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
306 void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayInt *n2oPtr)
308 static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
309 static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
310 if(!da || !mesh || !mfd)
311 throw MZCException("AppendMCFieldFrom : internal error !");
312 MEDFileFields *fs(mfd->getFields());
313 MEDFileMeshes *ms(mfd->getMeshes());
315 throw MZCException("AppendMCFieldFrom : internal error 2 !");
316 MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
317 DataArrayDouble *dadPtr(dad);
318 std::string fieldName;
321 fieldName=dadPtr->getName();
322 MCAuto<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
323 f->setName(fieldName);
328 MCAuto<DataArrayDouble> dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
332 MCAuto<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
333 MCAuto<MEDFileField1TS> f1ts(MEDFileField1TS::New());
334 f1ts->setFieldNoProfileSBT(f);
335 fmts->pushBackTimeStep(f1ts);
339 MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
340 DataArrayInt *daiPtr(dai);
343 fieldName=daiPtr->getName();
344 if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
346 MCAuto<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
347 f->setName(fieldName);
349 MCAuto<MEDFileIntFieldMultiTS> fmts(MEDFileIntFieldMultiTS::New());
350 MCAuto<MEDFileIntField1TS> f1ts(MEDFileIntField1TS::New());
352 f1ts->setFieldNoProfileSBT(f,daiPtr);
355 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
356 f1ts->setFieldNoProfileSBT(f,dai2);
358 fmts->pushBackTimeStep(f1ts);
362 else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
364 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
366 throw MZCException("AppendMCFieldFrom : internal error 3 !");
368 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daiPtr);
371 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
372 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
375 else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
377 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
379 throw MZCException("AppendMCFieldFrom : internal error 4 !");
381 mm->setFamilyFieldArr(1,daiPtr);
384 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
385 mm->setFamilyFieldArr(1,dai2);
391 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf)
394 throw MZCException("PutAtLevelDealOrder : internal error !");
395 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
396 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
398 throw MZCException("PutAtLevelDealOrder : internal error 2 !");
399 MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
400 mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
401 MCAuto<DataArrayInt> o2n(mesh->sortCellsInMEDFileFrmt());
402 const DataArrayInt *o2nPtr(o2n);
403 MCAuto<DataArrayInt> n2o;
404 mmu->setMeshAtLevel(meshDimRel,mesh);
405 const DataArrayInt *n2oPtr(0);
408 n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
410 if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
413 mm->setRenumFieldArr(meshDimRel,n2o);
416 std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
417 for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
419 MCAuto<DataArray> da(*it);
420 AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr);
424 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms)
427 throw MZCException("AssignSingleGTMeshes : internal error !");
428 MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
429 MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
431 throw MZCException("AssignSingleGTMeshes : internal error 2 !");
432 int meshDim(-std::numeric_limits<int>::max());
433 std::map<int, std::vector< MicroField > > ms2;
434 for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
436 const MEDCouplingUMesh *elt((*it).getMesh());
439 int myMeshDim(elt->getMeshDimension());
440 meshDim=std::max(meshDim,myMeshDim);
441 ms2[myMeshDim].push_back(*it);
446 for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
448 const std::vector< MicroField >& vs((*it).second);
451 PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0]);
455 MicroField merge(vs);
456 PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge);
461 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
464 throw MZCException("BuildCoordsFrom : internal error !");
465 vtkPoints *pts(ds->GetPoints());
467 throw MZCException("BuildCoordsFrom : internal error 2 !");
468 vtkDataArray *data(pts->GetData());
470 throw MZCException("BuildCoordsFrom : internal error 3 !");
471 MCAuto<DataArrayDouble> coords(ConvertVTKArrayToMCArrayDouble(data));
472 return coords.retn();
475 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa)
478 throw MZCException("AddNodeFields : internal error !");
479 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
480 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
482 throw MZCException("AddNodeFields : internal error 2 !");
483 MCAuto<MEDCouplingUMesh> mesh(mmu->getMeshAtLevel(0));
484 int nba(dsa->GetNumberOfArrays());
485 for(int i=0;i<nba;i++)
487 vtkDataArray *arr(dsa->GetArray(i));
488 const char *name(arr->GetName());
491 MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
493 AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,0);
497 std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayInt *part, vtkDataSetAttributes *dsa)
499 std::vector< MCAuto<DataArray> > ret;
502 int nba(dsa->GetNumberOfArrays());
503 for(int i=0;i<nba;i++)
505 vtkDataArray *arr(dsa->GetArray(i));
508 const char *name(arr->GetName());
509 int nbCompo(arr->GetNumberOfComponents());
510 vtkIdType nbTuples(arr->GetNumberOfTuples());
511 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
513 mcarr=mcarr->selectByTupleId(part->begin(),part->end());
514 mcarr->setName(name);
515 ret.push_back(mcarr);
520 std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
522 std::vector< MCAuto<DataArray> > ret;
525 int nba(dsa->GetNumberOfArrays());
526 for(int i=0;i<nba;i++)
528 vtkDataArray *arr(dsa->GetArray(i));
531 const char *name(arr->GetName());
532 int nbCompo(arr->GetNumberOfComponents());
533 vtkIdType nbTuples(arr->GetNumberOfTuples());
534 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
535 mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
536 mcarr->setName(name);
537 ret.push_back(mcarr);
542 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context)
545 throw MZCException("ConvertFromRectilinearGrid : internal error !");
547 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
548 ret->setMeshes(meshes);
549 MCAuto<MEDFileFields> fields(MEDFileFields::New());
550 ret->setFields(fields);
552 MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
553 meshes->pushMesh(cmesh);
554 MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
555 vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
558 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cx));
559 cmeshmc->setCoordsAt(0,arr);
563 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cy));
564 cmeshmc->setCoordsAt(1,arr);
568 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cz));
569 cmeshmc->setCoordsAt(2,arr);
571 std::string meshName(GetMeshNameWithContext(context));
572 cmeshmc->setName(meshName);
573 cmesh->setMesh(cmeshmc);
574 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
575 for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
577 MCAuto<DataArray> da(*it);
578 AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,0);
580 std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
581 for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
583 MCAuto<DataArray> da(*it);
584 AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,0);
588 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context)
591 throw MZCException("ConvertFromPolyData : internal error !");
593 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
594 ret->setMeshes(meshes);
595 MCAuto<MEDFileFields> fields(MEDFileFields::New());
596 ret->setFields(fields);
598 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
599 meshes->pushMesh(umesh);
600 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
601 umesh->setCoords(coords);
602 umesh->setName(GetMeshNameWithContext(context));
605 std::vector< MicroField > ms;
606 vtkCellArray *cd(ds->GetVerts());
609 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
610 if((const MEDCouplingUMesh *)subMesh)
612 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
613 offset+=subMesh->getNumberOfCells();
614 ms.push_back(MicroField(subMesh,cellFs));
617 vtkCellArray *cc(ds->GetLines());
620 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2));
621 if((const MEDCouplingUMesh *)subMesh)
623 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
624 offset+=subMesh->getNumberOfCells();
625 ms.push_back(MicroField(subMesh,cellFs));
628 vtkCellArray *cb(ds->GetPolys());
631 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
632 if((const MEDCouplingUMesh *)subMesh)
634 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
635 offset+=subMesh->getNumberOfCells();
636 ms.push_back(MicroField(subMesh,cellFs));
639 vtkCellArray *ca(ds->GetStrips());
642 MCAuto<DataArrayInt> ids;
643 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
644 if((const MEDCouplingUMesh *)subMesh)
646 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
647 offset+=subMesh->getNumberOfCells();
648 ms.push_back(MicroField(subMesh,cellFs));
651 AssignSingleGTMeshes(ret,ms);
652 AddNodeFields(ret,ds->GetPointData());
655 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context)
658 throw MZCException("ConvertFromUnstructuredGrid : internal error !");
660 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
661 ret->setMeshes(meshes);
662 MCAuto<MEDFileFields> fields(MEDFileFields::New());
663 ret->setFields(fields);
665 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
666 meshes->pushMesh(umesh);
667 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
668 umesh->setCoords(coords);
669 umesh->setName(GetMeshNameWithContext(context));
670 vtkIdType nbCells(ds->GetNumberOfCells());
671 vtkCellArray *ca(ds->GetCells());
674 vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries());
675 vtkIdType *caPtr(ca->GetPointer());
676 vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
678 throw MZCException("ConvertFromUnstructuredGrid : internal error");
679 vtkIdTypeArray *cla(ds->GetCellLocationsArray());
680 const vtkIdType *claPtr(cla->GetPointer(0));
682 throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
683 const unsigned char *ctPtr(ct->GetPointer(0));
684 std::map<int,int> m(ComputeMapOfType());
685 MCAuto<DataArrayInt> lev(DataArrayInt::New()) ; lev->alloc(nbCells,1);
686 int *levPtr(lev->getPointer());
687 for(vtkIdType i=0;i<nbCells;i++)
689 std::map<int,int>::iterator it(m.find(ctPtr[i]));
692 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
693 levPtr[i]=cm.getDimension();
697 std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
698 throw MZCException(oss.str());
702 int meshDim(lev->getMaxValue(dummy));
703 MCAuto<DataArrayInt> levs(lev->getDifferentValues());
704 std::vector< MicroField > ms;
705 vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations());
706 for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
708 MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
709 m0->setCoords(coords); m0->allocateCells();
710 MCAuto<DataArrayInt> cellIdsCurLev(lev->findIdsEqual(*curLev));
711 for(const int *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
713 std::map<int,int>::iterator it(m.find(ctPtr[*cellId]));
714 vtkIdType offset(claPtr[*cellId]);
715 vtkIdType sz(caPtr[offset]);
716 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)(*it).second);
717 if(ct!=INTERP_KERNEL::NORM_POLYHED)
719 std::vector<int> conn2(sz);
720 for(int kk=0;kk<sz;kk++)
721 conn2[kk]=caPtr[offset+1+kk];
722 m0->insertNextCell(ct,sz,&conn2[0]);
726 if(!faces || !faceLoc)
727 throw MZCException("ConvertFromUnstructuredGrid : faces are expected when there are polyhedra !");
728 const vtkIdType *facPtr(faces->GetPointer(0)),*facLocPtr(faceLoc->GetPointer(0));
729 std::vector<int> conn;
730 int off0(facLocPtr[*cellId]);
731 int nbOfFaces(facPtr[off0++]);
732 for(int k=0;k<nbOfFaces;k++)
734 int nbOfNodesInFace(facPtr[off0++]);
735 std::copy(facPtr+off0,facPtr+off0+nbOfNodesInFace,std::back_inserter(conn));
736 off0+=nbOfNodesInFace;
740 m0->insertNextCell(ct,conn.size(),&conn[0]);
743 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
744 ms.push_back(MicroField(m0,cellFs));
746 AssignSingleGTMeshes(ret,ms);
747 AddNodeFields(ret,ds->GetPointData());
752 vtkMEDWriter::vtkMEDWriter():WriteAllTimeSteps(0),FileName(0),IsTouched(false)
756 vtkMEDWriter::~vtkMEDWriter()
760 vtkInformationDataObjectMetaDataKey *GetMEDReaderMetaDataIfAny()
762 static const char ZE_KEY[]="vtkMEDReader::META_DATA";
763 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
764 if(!gd->hasKey(ZE_KEY))
766 std::string ptSt(gd->value(ZE_KEY));
768 std::istringstream iss(ptSt); iss >> pt;
769 return reinterpret_cast<vtkInformationDataObjectMetaDataKey *>(pt);
772 bool IsInformationOK(vtkInformation *info)
774 vtkInformationDataObjectMetaDataKey *key(GetMEDReaderMetaDataIfAny());
777 // Check the information contain meta data key
781 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
785 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
786 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
789 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
791 vtkStdString &st(verticesNames2->GetValue(i));
792 if(st=="MeshesFamsGrps")
801 Grp(const std::string& name):_name(name) { }
802 void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
803 std::string getName() const { return _name; }
804 std::vector<std::string> getFamilies() const { return _fams; }
807 std::vector<std::string> _fams;
813 Fam(const std::string& name)
815 static const char ZE_SEP[]="@@][@@";
816 std::size_t pos(name.find(ZE_SEP));
817 std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
818 std::istringstream iss(name1);
822 std::string getName() const { return _name; }
823 int getID() const { return _id; }
829 void LoadFamGrpMapInfo(vtkMutableDirectedGraph *sil, std::string& meshName, std::vector<Grp>& groups, std::vector<Fam>& fams)
832 throw MZCException("LoadFamGrpMapInfo : internal error !");
834 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
835 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
838 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
840 vtkStdString &st(verticesNames2->GetValue(i));
841 if(st=="MeshesFamsGrps")
848 throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
849 vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
850 sil->GetAdjacentVertices(id0,it0);
852 while(it0->HasNext())
854 vtkIdType id1(it0->Next());
855 std::string mName((const char *)verticesNames2->GetValue(id1));
857 vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
858 sil->GetAdjacentVertices(id1,it1);
859 vtkIdType idZeGrps(it1->Next());//zeGroups
860 vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
861 sil->GetAdjacentVertices(idZeGrps,itGrps);
862 while(itGrps->HasNext())
864 vtkIdType idg(itGrps->Next());
865 Grp grp((const char *)verticesNames2->GetValue(idg));
866 vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
867 sil->GetAdjacentVertices(idg,itGrps2);
868 std::vector<std::string> famsOnGroup;
869 while(itGrps2->HasNext())
871 vtkIdType idgf(itGrps2->Next());
872 famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
874 grp.setFamilies(famsOnGroup);
876 groups.push_back(grp);
879 vtkIdType idZeFams(it1->Next());//zeFams
881 vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
882 sil->GetAdjacentVertices(idZeFams,itFams);
883 while(itFams->HasNext())
885 vtkIdType idf(itFams->Next());
886 Fam fam((const char *)verticesNames2->GetValue(idf));
894 int vtkMEDWriter::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
896 //std::cerr << "########################################## vtkMEDWriter::RequestInformation ########################################## " << (const char *) this->FileName << std::endl;
900 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context)
903 throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
904 vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
905 vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
906 vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
909 ConvertFromPolyData(mfd,ds2,context);
913 ConvertFromUnstructuredGrid(mfd,ds3,context);
917 ConvertFromRectilinearGrid(mfd,ds4,context);
920 throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
923 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context)
926 throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
927 int nbBlocks(ds->GetNumberOfBlocks());
928 if(nbBlocks==1 && context.empty())
930 vtkDataObject *uniqueElt(ds->GetBlock(0));
932 throw MZCException("Unique elt in multiblock is NULL !");
933 vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
936 WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context);
940 for(int i=0;i<nbBlocks;i++)
942 vtkDataObject *elt(ds->GetBlock(i));
943 std::vector<int> context2;
944 context2.push_back(i);
947 std::ostringstream oss; oss << "In context ";
948 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
949 oss << " at pos #" << i << " elt is NULL !";
950 throw MZCException(oss.str());
952 vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
955 WriteMEDFileFromVTKDataSet(mfd,elt1,context);
958 vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
961 WriteMEDFileFromVTKMultiBlock(mfd,elt2,context);
964 std::ostringstream oss; oss << "In context ";
965 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
966 oss << " at pos #" << i << " elt not recognized data type !";
967 throw MZCException(oss.str());
971 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input)
974 throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
975 std::vector<int> context;
976 vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
979 WriteMEDFileFromVTKDataSet(mfd,input1,context);
982 vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
985 WriteMEDFileFromVTKMultiBlock(mfd,input2,context);
988 throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
991 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
997 MEDFileMeshes *meshes(mfd->getMeshes());
1000 if(meshes->getNumberOfMeshes()!=1)
1002 MEDFileMesh *mm(meshes->getMeshAtPos(0));
1005 mm->setName(meshName);
1006 for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
1007 mm->setFamilyId((*it).getName(),(*it).getID());
1008 for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
1009 mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
1010 MEDFileFields *fields(mfd->getFields());
1013 for(int i=0;i<fields->getNumberOfFields();i++)
1015 MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
1018 fmts->setMeshName(meshName);
1022 int vtkMEDWriter::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
1024 //std::cerr << "########################################## vtkMEDWriter::RequestData ########################################## " << (const char *) this->FileName << std::endl;
1027 vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
1028 std::string meshName;
1029 std::vector<Grp> groups; std::vector<Fam> fams;
1030 if(IsInformationOK(inputInfo))
1032 vtkMutableDirectedGraph *famGrpGraph(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
1033 LoadFamGrpMapInfo(famGrpGraph,meshName,groups,fams);
1035 vtkInformation *outInfo(outputVector->GetInformationObject(0));
1036 vtkDataObject *input(vtkDataObject::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
1038 throw MZCException("Not recognized data object in input of the MEDWriter ! Maybe not implemented yet !");
1039 MCAuto<MEDFileData> mfd(MEDFileData::New());
1040 WriteMEDFileFromVTKGDS(mfd,input);
1041 PutFamGrpInfoIfAny(mfd,meshName,groups,fams);
1042 mfd->write(this->FileName,this->IsTouched?0:2); this->IsTouched=true;
1043 outInfo->Set(vtkDataObject::DATA_OBJECT(),input);
1045 catch(MZCException& e)
1047 std::ostringstream oss;
1048 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;
1049 if(this->HasObserver("ErrorEvent") )
1050 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
1052 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
1053 vtkObject::BreakOnError();
1059 void vtkMEDWriter::PrintSelf(ostream& os, vtkIndent indent)
1061 this->Superclass::PrintSelf(os, indent);
1064 int vtkMEDWriter::Write()