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 "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"
60 #include "MEDFileMesh.hxx"
61 #include "MEDFileField.hxx"
62 #include "MEDFileData.hxx"
63 #include "MEDCouplingMemArray.hxx"
64 #include "MEDCouplingFieldInt.hxx"
65 #include "MEDCouplingFieldDouble.hxx"
66 #include "MEDCouplingRefCountObject.hxx"
72 using MEDCoupling::MEDFileData;
73 using MEDCoupling::MEDFileMesh;
74 using MEDCoupling::MEDFileCMesh;
75 using MEDCoupling::MEDFileUMesh;
76 using MEDCoupling::MEDFileFields;
77 using MEDCoupling::MEDFileMeshes;
79 using MEDCoupling::MEDFileIntField1TS;
80 using MEDCoupling::MEDFileField1TS;
81 using MEDCoupling::MEDFileIntFieldMultiTS;
82 using MEDCoupling::MEDFileFieldMultiTS;
83 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
84 using MEDCoupling::DataArray;
85 using MEDCoupling::DataArrayInt;
86 using MEDCoupling::DataArrayDouble;
87 using MEDCoupling::MEDCouplingMesh;
88 using MEDCoupling::MEDCouplingUMesh;
89 using MEDCoupling::MEDCouplingCMesh;
90 using MEDCoupling::MEDCouplingFieldDouble;
91 using MEDCoupling::MEDCouplingFieldInt;
92 using MEDCoupling::MCAuto;
94 vtkStandardNewMacro(vtkMEDWriter);
98 class MZCException : public std::exception
101 MZCException(const std::string& s):_reason(s) { }
102 virtual const char *what() const throw() { return _reason.c_str(); }
103 virtual ~MZCException() throw() { }
110 std::map<int,int> ComputeMapOfType()
112 std::map<int,int> ret;
113 int nbOfTypesInMC(sizeof(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER)/sizeof(int));
114 for(int i=0;i<nbOfTypesInMC;i++)
116 int vtkId(MEDCouplingUMesh::MEDCOUPLING2VTKTYPETRADUCER[i]);
123 std::string GetMeshNameWithContext(const std::vector<int>& context)
125 static const char DFT_MESH_NAME[]="Mesh";
127 return DFT_MESH_NAME;
128 std::ostringstream oss; oss << DFT_MESH_NAME;
129 for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
134 DataArrayInt *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
137 throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
138 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
139 std::size_t nbElts(nbTuples*nbComp);
140 MCAuto<DataArrayInt> ret(DataArrayInt::New());
141 ret->alloc(nbTuples,nbComp);
142 for(int i=0;i<nbComp;i++)
144 const char *comp(data->GetComponentName(i));
146 ret->setInfoOnComponent(i,comp);
148 int *ptOut(ret->getPointer());
149 vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
152 const int *pt(d0->GetPointer(0));
153 std::copy(pt,pt+nbElts,ptOut);
156 vtkLongArray *d1(vtkLongArray::SafeDownCast(data));
159 const long *pt(d1->GetPointer(0));
160 std::copy(pt,pt+nbElts,ptOut);
163 std::ostringstream oss;
164 oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
165 throw MZCException(oss.str());
168 DataArrayDouble *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
171 throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
172 int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
173 std::size_t nbElts(nbTuples*nbComp);
174 MCAuto<DataArrayDouble> ret(DataArrayDouble::New());
175 ret->alloc(nbTuples,nbComp);
176 for(int i=0;i<nbComp;i++)
178 const char *comp(data->GetComponentName(i));
180 ret->setInfoOnComponent(i,comp);
182 double *ptOut(ret->getPointer());
183 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
186 const float *pt(d0->GetPointer(0));
187 for(std::size_t i=0;i<nbElts;i++)
191 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
194 const double *pt(d1->GetPointer(0));
195 std::copy(pt,pt+nbElts,ptOut);
198 std::ostringstream oss;
199 oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << typeid(*data).name() << "\" type !";
200 throw MZCException(oss.str());
203 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
206 throw MZCException("ConvertVTKArrayToMCArray : internal error !");
207 vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
208 vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
210 return ConvertVTKArrayToMCArrayDouble(data);
211 vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
212 vtkLongArray *d3(vtkLongArray::SafeDownCast(data));
214 return ConvertVTKArrayToMCArrayInt(data);
215 std::ostringstream oss;
216 oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
217 throw MZCException(oss.str());
220 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
222 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
223 subMesh->setCoords(coords); subMesh->allocateCells();
224 int nbCells(ca->GetNumberOfCells());
227 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
228 const vtkIdType *conn(ca->GetPointer());
229 for(int i=0;i<nbCells;i++)
232 std::vector<int> conn2(sz);
233 for(int jj=0;jj<sz;jj++)
235 subMesh->insertNextCell(type,sz,&conn2[0]);
238 return subMesh.retn();
241 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayInt>& ids)
243 MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
244 subMesh->setCoords(coords); subMesh->allocateCells();
245 int nbCells(ca->GetNumberOfCells());
248 vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
249 const vtkIdType *conn(ca->GetPointer());
250 ids=DataArrayInt::New() ; ids->alloc(0,1);
251 for(int i=0;i<nbCells;i++)
257 for(int j=0;j<nbTri;j++,conn++)
259 int conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
260 subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
261 ids->pushBackSilent(i);
266 std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
267 throw MZCException(oss.str());
271 return subMesh.retn();
277 MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
278 MicroField(const std::vector< MicroField >& vs);
279 void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
280 MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
281 std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
283 MCAuto<MEDCouplingUMesh> _m;
284 std::vector<MCAuto<DataArray> > _cellFs;
285 std::vector<MCAuto<DataArray> > _nodeFs;
288 MicroField::MicroField(const std::vector< MicroField >& vs)
290 std::size_t sz(vs.size());
291 std::vector<const MEDCouplingUMesh *> vs2(sz);
292 std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
294 for(std::size_t ii=0;ii<sz;ii++)
296 vs2[ii]=vs[ii].getMesh();
297 arrs2[ii]=vs[ii].getCellFields();
299 nbElts=arrs2[ii].size();
301 if(arrs2[ii].size()!=nbElts)
302 throw MZCException("MicroField cstor : internal error !");
304 _cellFs.resize(nbElts);
305 for(int ii=0;ii<nbElts;ii++)
307 std::vector<const DataArray *> arrsTmp(sz);
308 for(std::size_t jj=0;jj<sz;jj++)
310 arrsTmp[jj]=arrs2[jj][ii];
312 _cellFs[ii]=DataArray::Aggregate(arrsTmp);
314 _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
317 void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayInt *n2oPtr)
319 static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
320 static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
321 if(!da || !mesh || !mfd)
322 throw MZCException("AppendMCFieldFrom : internal error !");
323 MEDFileFields *fs(mfd->getFields());
324 MEDFileMeshes *ms(mfd->getMeshes());
326 throw MZCException("AppendMCFieldFrom : internal error 2 !");
327 MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
328 DataArrayDouble *dadPtr(dad);
329 std::string fieldName;
332 fieldName=dadPtr->getName();
333 MCAuto<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
334 f->setName(fieldName);
339 MCAuto<DataArrayDouble> dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
343 MCAuto<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
344 MCAuto<MEDFileField1TS> f1ts(MEDFileField1TS::New());
345 f1ts->setFieldNoProfileSBT(f);
346 fmts->pushBackTimeStep(f1ts);
350 MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
351 DataArrayInt *daiPtr(dai);
354 fieldName=daiPtr->getName();
355 if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
357 MCAuto<MEDCouplingFieldInt> f(MEDCouplingFieldInt::New(tf));
358 f->setName(fieldName);
360 MCAuto<MEDFileIntFieldMultiTS> fmts(MEDFileIntFieldMultiTS::New());
361 MCAuto<MEDFileIntField1TS> f1ts(MEDFileIntField1TS::New());
365 f1ts->setFieldNoProfileSBT(f);
369 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
371 f1ts->setFieldNoProfileSBT(f);
373 fmts->pushBackTimeStep(f1ts);
377 else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
379 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
381 throw MZCException("AppendMCFieldFrom : internal error 3 !");
383 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daiPtr);
386 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
387 mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
390 else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
392 MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
394 throw MZCException("AppendMCFieldFrom : internal error 4 !");
396 mm->setFamilyFieldArr(1,daiPtr);
399 MCAuto<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
400 mm->setFamilyFieldArr(1,dai2);
406 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf)
409 throw MZCException("PutAtLevelDealOrder : internal error !");
410 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
411 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
413 throw MZCException("PutAtLevelDealOrder : internal error 2 !");
414 MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
415 mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
416 MCAuto<DataArrayInt> o2n(mesh->sortCellsInMEDFileFrmt());
417 const DataArrayInt *o2nPtr(o2n);
418 MCAuto<DataArrayInt> n2o;
419 mmu->setMeshAtLevel(meshDimRel,mesh);
420 const DataArrayInt *n2oPtr(0);
423 n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
425 if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
428 mm->setRenumFieldArr(meshDimRel,n2o);
431 std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
432 for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
434 MCAuto<DataArray> da(*it);
435 AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr);
439 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms)
442 throw MZCException("AssignSingleGTMeshes : internal error !");
443 MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
444 MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
446 throw MZCException("AssignSingleGTMeshes : internal error 2 !");
447 int meshDim(-std::numeric_limits<int>::max());
448 std::map<int, std::vector< MicroField > > ms2;
449 for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
451 const MEDCouplingUMesh *elt((*it).getMesh());
454 int myMeshDim(elt->getMeshDimension());
455 meshDim=std::max(meshDim,myMeshDim);
456 ms2[myMeshDim].push_back(*it);
461 for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
463 const std::vector< MicroField >& vs((*it).second);
466 PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0]);
470 MicroField merge(vs);
471 PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge);
476 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
479 throw MZCException("BuildCoordsFrom : internal error !");
480 vtkPoints *pts(ds->GetPoints());
482 throw MZCException("BuildCoordsFrom : internal error 2 !");
483 vtkDataArray *data(pts->GetData());
485 throw MZCException("BuildCoordsFrom : internal error 3 !");
486 MCAuto<DataArrayDouble> coords(ConvertVTKArrayToMCArrayDouble(data));
487 return coords.retn();
490 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa)
493 throw MZCException("AddNodeFields : internal error !");
494 MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
495 MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
497 throw MZCException("AddNodeFields : internal error 2 !");
498 MCAuto<MEDCouplingUMesh> mesh;
499 if(!mmu->getNonEmptyLevels().empty())
500 mesh=mmu->getMeshAtLevel(0);
503 mesh=MEDCouplingUMesh::Build0DMeshFromCoords(mmu->getCoords());
504 mesh->setName(mmu->getName());
506 int nba(dsa->GetNumberOfArrays());
507 for(int i=0;i<nba;i++)
509 vtkDataArray *arr(dsa->GetArray(i));
510 const char *name(arr->GetName());
513 MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
515 AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,0);
519 std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayInt *part, vtkDataSetAttributes *dsa)
521 std::vector< MCAuto<DataArray> > ret;
524 int nba(dsa->GetNumberOfArrays());
525 for(int i=0;i<nba;i++)
527 vtkDataArray *arr(dsa->GetArray(i));
530 const char *name(arr->GetName());
531 int nbCompo(arr->GetNumberOfComponents());
532 vtkIdType nbTuples(arr->GetNumberOfTuples());
533 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
535 mcarr=mcarr->selectByTupleId(part->begin(),part->end());
536 mcarr->setName(name);
537 ret.push_back(mcarr);
542 std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
544 std::vector< MCAuto<DataArray> > ret;
547 int nba(dsa->GetNumberOfArrays());
548 for(int i=0;i<nba;i++)
550 vtkDataArray *arr(dsa->GetArray(i));
553 const char *name(arr->GetName());
554 int nbCompo(arr->GetNumberOfComponents());
555 vtkIdType nbTuples(arr->GetNumberOfTuples());
556 MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
557 mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
558 mcarr->setName(name);
559 ret.push_back(mcarr);
564 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context)
567 throw MZCException("ConvertFromRectilinearGrid : internal error !");
569 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
570 ret->setMeshes(meshes);
571 MCAuto<MEDFileFields> fields(MEDFileFields::New());
572 ret->setFields(fields);
574 MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
575 meshes->pushMesh(cmesh);
576 MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
577 vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
580 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cx));
581 cmeshmc->setCoordsAt(0,arr);
585 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cy));
586 cmeshmc->setCoordsAt(1,arr);
590 MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cz));
591 cmeshmc->setCoordsAt(2,arr);
593 std::string meshName(GetMeshNameWithContext(context));
594 cmeshmc->setName(meshName);
595 cmesh->setMesh(cmeshmc);
596 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
597 for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
599 MCAuto<DataArray> da(*it);
600 AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,0);
602 std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
603 for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
605 MCAuto<DataArray> da(*it);
606 AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,0);
610 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context)
613 throw MZCException("ConvertFromPolyData : internal error !");
615 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
616 ret->setMeshes(meshes);
617 MCAuto<MEDFileFields> fields(MEDFileFields::New());
618 ret->setFields(fields);
620 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
621 meshes->pushMesh(umesh);
622 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
623 umesh->setCoords(coords);
624 umesh->setName(GetMeshNameWithContext(context));
627 std::vector< MicroField > ms;
628 vtkCellArray *cd(ds->GetVerts());
631 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
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 *cc(ds->GetLines());
642 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2));
643 if((const MEDCouplingUMesh *)subMesh)
645 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
646 offset+=subMesh->getNumberOfCells();
647 ms.push_back(MicroField(subMesh,cellFs));
650 vtkCellArray *cb(ds->GetPolys());
653 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
654 if((const MEDCouplingUMesh *)subMesh)
656 std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
657 offset+=subMesh->getNumberOfCells();
658 ms.push_back(MicroField(subMesh,cellFs));
661 vtkCellArray *ca(ds->GetStrips());
664 MCAuto<DataArrayInt> ids;
665 MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
666 if((const MEDCouplingUMesh *)subMesh)
668 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
669 offset+=subMesh->getNumberOfCells();
670 ms.push_back(MicroField(subMesh,cellFs));
673 AssignSingleGTMeshes(ret,ms);
674 AddNodeFields(ret,ds->GetPointData());
677 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context)
680 throw MZCException("ConvertFromUnstructuredGrid : internal error !");
682 MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
683 ret->setMeshes(meshes);
684 MCAuto<MEDFileFields> fields(MEDFileFields::New());
685 ret->setFields(fields);
687 MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
688 meshes->pushMesh(umesh);
689 MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
690 umesh->setCoords(coords);
691 umesh->setName(GetMeshNameWithContext(context));
692 vtkIdType nbCells(ds->GetNumberOfCells());
693 vtkCellArray *ca(ds->GetCells());
696 vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries());
697 vtkIdType *caPtr(ca->GetPointer());
698 vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
700 throw MZCException("ConvertFromUnstructuredGrid : internal error");
701 vtkIdTypeArray *cla(ds->GetCellLocationsArray());
702 const vtkIdType *claPtr(cla->GetPointer(0));
704 throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
705 const unsigned char *ctPtr(ct->GetPointer(0));
706 std::map<int,int> m(ComputeMapOfType());
707 MCAuto<DataArrayInt> lev(DataArrayInt::New()) ; lev->alloc(nbCells,1);
708 int *levPtr(lev->getPointer());
709 for(vtkIdType i=0;i<nbCells;i++)
711 std::map<int,int>::iterator it(m.find(ctPtr[i]));
714 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
715 levPtr[i]=cm.getDimension();
719 std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
720 throw MZCException(oss.str());
724 MCAuto<DataArrayInt> levs(lev->getDifferentValues());
725 std::vector< MicroField > ms;
726 vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations());
727 for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
729 MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
730 m0->setCoords(coords); m0->allocateCells();
731 MCAuto<DataArrayInt> cellIdsCurLev(lev->findIdsEqual(*curLev));
732 for(const int *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
734 std::map<int,int>::iterator it(m.find(ctPtr[*cellId]));
735 vtkIdType offset(claPtr[*cellId]);
736 vtkIdType sz(caPtr[offset]);
737 INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)(*it).second);
738 if(ct!=INTERP_KERNEL::NORM_POLYHED)
740 std::vector<int> conn2(sz);
741 for(int kk=0;kk<sz;kk++)
742 conn2[kk]=caPtr[offset+1+kk];
743 m0->insertNextCell(ct,sz,&conn2[0]);
747 if(!faces || !faceLoc)
748 throw MZCException("ConvertFromUnstructuredGrid : faces are expected when there are polyhedra !");
749 const vtkIdType *facPtr(faces->GetPointer(0)),*facLocPtr(faceLoc->GetPointer(0));
750 std::vector<int> conn;
751 int off0(facLocPtr[*cellId]);
752 int nbOfFaces(facPtr[off0++]);
753 for(int k=0;k<nbOfFaces;k++)
755 int nbOfNodesInFace(facPtr[off0++]);
756 std::copy(facPtr+off0,facPtr+off0+nbOfNodesInFace,std::back_inserter(conn));
757 off0+=nbOfNodesInFace;
761 m0->insertNextCell(ct,conn.size(),&conn[0]);
764 std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
765 ms.push_back(MicroField(m0,cellFs));
767 AssignSingleGTMeshes(ret,ms);
768 AddNodeFields(ret,ds->GetPointData());
773 vtkMEDWriter::vtkMEDWriter():WriteAllTimeSteps(0),FileName(0),IsTouched(false)
777 vtkMEDWriter::~vtkMEDWriter()
781 vtkInformationDataObjectMetaDataKey *GetMEDReaderMetaDataIfAny()
783 static const char ZE_KEY[]="vtkMEDReader::META_DATA";
784 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
785 if(!gd->hasKey(ZE_KEY))
787 std::string ptSt(gd->value(ZE_KEY));
789 std::istringstream iss(ptSt); iss >> pt;
790 return reinterpret_cast<vtkInformationDataObjectMetaDataKey *>(pt);
793 bool IsInformationOK(vtkInformation *info)
795 vtkInformationDataObjectMetaDataKey *key(GetMEDReaderMetaDataIfAny());
798 // Check the information contain meta data key
802 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
806 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
807 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
810 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
812 vtkStdString &st(verticesNames2->GetValue(i));
813 if(st=="MeshesFamsGrps")
822 Grp(const std::string& name):_name(name) { }
823 void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
824 std::string getName() const { return _name; }
825 std::vector<std::string> getFamilies() const { return _fams; }
828 std::vector<std::string> _fams;
834 Fam(const std::string& name)
836 static const char ZE_SEP[]="@@][@@";
837 std::size_t pos(name.find(ZE_SEP));
838 std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
839 std::istringstream iss(name1);
843 std::string getName() const { return _name; }
844 int getID() const { return _id; }
850 void LoadFamGrpMapInfo(vtkMutableDirectedGraph *sil, std::string& meshName, std::vector<Grp>& groups, std::vector<Fam>& fams)
853 throw MZCException("LoadFamGrpMapInfo : internal error !");
855 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
856 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
859 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
861 vtkStdString &st(verticesNames2->GetValue(i));
862 if(st=="MeshesFamsGrps")
869 throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
870 vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
871 sil->GetAdjacentVertices(id0,it0);
873 while(it0->HasNext())
875 vtkIdType id1(it0->Next());
876 std::string mName((const char *)verticesNames2->GetValue(id1));
878 vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
879 sil->GetAdjacentVertices(id1,it1);
880 vtkIdType idZeGrps(it1->Next());//zeGroups
881 vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
882 sil->GetAdjacentVertices(idZeGrps,itGrps);
883 while(itGrps->HasNext())
885 vtkIdType idg(itGrps->Next());
886 Grp grp((const char *)verticesNames2->GetValue(idg));
887 vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
888 sil->GetAdjacentVertices(idg,itGrps2);
889 std::vector<std::string> famsOnGroup;
890 while(itGrps2->HasNext())
892 vtkIdType idgf(itGrps2->Next());
893 famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
895 grp.setFamilies(famsOnGroup);
897 groups.push_back(grp);
900 vtkIdType idZeFams(it1->Next());//zeFams
902 vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
903 sil->GetAdjacentVertices(idZeFams,itFams);
904 while(itFams->HasNext())
906 vtkIdType idf(itFams->Next());
907 Fam fam((const char *)verticesNames2->GetValue(idf));
915 int vtkMEDWriter::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
917 //std::cerr << "########################################## vtkMEDWriter::RequestInformation ########################################## " << (const char *) this->FileName << std::endl;
921 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context)
924 throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
925 vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
926 vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
927 vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
930 ConvertFromPolyData(mfd,ds2,context);
934 ConvertFromUnstructuredGrid(mfd,ds3,context);
938 ConvertFromRectilinearGrid(mfd,ds4,context);
941 throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
944 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context)
947 throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
948 int nbBlocks(ds->GetNumberOfBlocks());
949 if(nbBlocks==1 && context.empty())
951 vtkDataObject *uniqueElt(ds->GetBlock(0));
953 throw MZCException("Unique elt in multiblock is NULL !");
954 vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
957 WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context);
961 for(int i=0;i<nbBlocks;i++)
963 vtkDataObject *elt(ds->GetBlock(i));
964 std::vector<int> context2;
965 context2.push_back(i);
968 std::ostringstream oss; oss << "In context ";
969 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
970 oss << " at pos #" << i << " elt is NULL !";
971 throw MZCException(oss.str());
973 vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
976 WriteMEDFileFromVTKDataSet(mfd,elt1,context);
979 vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
982 WriteMEDFileFromVTKMultiBlock(mfd,elt2,context);
985 std::ostringstream oss; oss << "In context ";
986 std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
987 oss << " at pos #" << i << " elt not recognized data type !";
988 throw MZCException(oss.str());
992 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input)
995 throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
996 std::vector<int> context;
997 vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
1000 WriteMEDFileFromVTKDataSet(mfd,input1,context);
1003 vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
1006 WriteMEDFileFromVTKMultiBlock(mfd,input2,context);
1009 throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
1012 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
1016 if(meshName.empty())
1018 MEDFileMeshes *meshes(mfd->getMeshes());
1021 if(meshes->getNumberOfMeshes()!=1)
1023 MEDFileMesh *mm(meshes->getMeshAtPos(0));
1026 mm->setName(meshName);
1027 for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
1028 mm->setFamilyId((*it).getName(),(*it).getID());
1029 for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
1030 mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
1031 MEDFileFields *fields(mfd->getFields());
1034 for(int i=0;i<fields->getNumberOfFields();i++)
1036 MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
1039 fmts->setMeshName(meshName);
1043 int vtkMEDWriter::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
1045 //std::cerr << "########################################## vtkMEDWriter::RequestData ########################################## " << (const char *) this->FileName << std::endl;
1048 vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
1049 std::string meshName;
1050 std::vector<Grp> groups; std::vector<Fam> fams;
1051 if(IsInformationOK(inputInfo))
1053 vtkMutableDirectedGraph *famGrpGraph(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
1054 LoadFamGrpMapInfo(famGrpGraph,meshName,groups,fams);
1056 vtkInformation *outInfo(outputVector->GetInformationObject(0));
1057 vtkDataObject *input(vtkDataObject::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
1059 throw MZCException("Not recognized data object in input of the MEDWriter ! Maybe not implemented yet !");
1060 MCAuto<MEDFileData> mfd(MEDFileData::New());
1061 WriteMEDFileFromVTKGDS(mfd,input);
1062 PutFamGrpInfoIfAny(mfd,meshName,groups,fams);
1063 mfd->write(this->FileName,this->IsTouched?0:2); this->IsTouched=true;
1064 outInfo->Set(vtkDataObject::DATA_OBJECT(),input);
1066 catch(MZCException& e)
1068 std::ostringstream oss;
1069 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;
1070 if(this->HasObserver("ErrorEvent") )
1071 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
1073 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
1074 vtkObject::BreakOnError();
1080 void vtkMEDWriter::PrintSelf(ostream& os, vtkIndent indent)
1082 this->Superclass::PrintSelf(os, indent);
1085 int vtkMEDWriter::Write()