]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/MEDWriter/IO/vtkMEDWriter.cxx
Salome HOME
MEDWriter plugin is here.
[modules/paravis.git] / src / Plugins / MEDWriter / IO / vtkMEDWriter.cxx
1 // Copyright (C) 2016  CEA/DEN, EDF R&D
2 //
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.
7 //
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.
12 //
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
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 // Author : Anthony Geay (EDF R&D)
20
21 #include "vtkMEDWriter.h"
22
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkIntArray.h"
25 #include "vtkCellData.h"
26 #include "vtkPointData.h"
27 #include "vtkFloatArray.h"
28
29 #include "vtkStreamingDemandDrivenPipeline.h"
30 #include "vtkInformationDataObjectMetaDataKey.h"
31 #include "vtkUnstructuredGrid.h"
32 #include "vtkMultiBlockDataSet.h"
33 #include "vtkRectilinearGrid.h"
34 #include "vtkInformationStringKey.h"
35 #include "vtkAlgorithmOutput.h"
36 #include "vtkObjectFactory.h"
37 #include "vtkMutableDirectedGraph.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPolyData.h"
40 #include "vtkDataSet.h"
41 #include "vtkInformationVector.h"
42 #include "vtkInformation.h"
43 #include "vtkDataArraySelection.h"
44 #include "vtkTimeStamp.h"
45 #include "vtkInEdgeIterator.h"
46 #include "vtkInformationDataObjectKey.h"
47 #include "vtkExecutive.h"
48 #include "vtkVariantArray.h"
49 #include "vtkStringArray.h"
50 #include "vtkDoubleArray.h"
51 #include "vtkCharArray.h"
52 #include "vtkUnsignedCharArray.h"
53 #include "vtkDataSetAttributes.h"
54 #include "vtkDemandDrivenPipeline.h"
55 #include "vtkDataObjectTreeIterator.h"
56 #include "vtkWarpScalar.h"
57
58 #include "MEDFileMesh.hxx"
59 #include "MEDFileField.hxx"
60 #include "MEDFileData.hxx"
61 #include "MEDCouplingMemArray.hxx"
62 #include "MEDCouplingFieldDouble.hxx"
63 #include "MEDCouplingAutoRefCountObjectPtr.hxx"
64
65 #include <map>
66 #include <deque>
67 #include <sstream>
68
69 using ParaMEDMEM::MEDFileData;
70 using ParaMEDMEM::MEDFileMesh;
71 using ParaMEDMEM::MEDFileCMesh;
72 using ParaMEDMEM::MEDFileUMesh;
73 using ParaMEDMEM::MEDFileFields;
74 using ParaMEDMEM::MEDFileMeshes;
75
76 using ParaMEDMEM::MEDFileIntField1TS;
77 using ParaMEDMEM::MEDFileField1TS;
78 using ParaMEDMEM::MEDFileIntFieldMultiTS;
79 using ParaMEDMEM::MEDFileFieldMultiTS;
80 using ParaMEDMEM::MEDFileAnyTypeFieldMultiTS;
81 using ParaMEDMEM::DataArray;
82 using ParaMEDMEM::DataArrayInt;
83 using ParaMEDMEM::DataArrayDouble;
84 using ParaMEDMEM::MEDCouplingMesh;
85 using ParaMEDMEM::MEDCouplingUMesh;
86 using ParaMEDMEM::MEDCouplingCMesh;
87 using ParaMEDMEM::MEDCouplingFieldDouble;
88 using ParaMEDMEM::MEDCouplingAutoRefCountObjectPtr;
89
90 vtkStandardNewMacro(vtkMEDWriter);
91
92 ///////////////////
93
94 class MZCException : public std::exception
95 {
96 public:
97   MZCException(const std::string& s):_reason(s) { }
98   virtual const char *what() const throw() { return _reason.c_str(); }
99   virtual ~MZCException() throw() { }
100 private:
101   std::string _reason;
102 };
103
104 ///////////////////
105
106 std::map<int,int> ComputeMapOfType()
107 {
108   std::map<int,int> ret;
109   int nbOfTypesInMC(sizeof(MEDCouplingUMesh::PARAMEDMEM2VTKTYPETRADUCER)/sizeof(int));
110   for(int i=0;i<nbOfTypesInMC;i++)
111     {
112       int vtkId(MEDCouplingUMesh::PARAMEDMEM2VTKTYPETRADUCER[i]);
113       if(vtkId!=-1)
114         ret[vtkId]=i;
115     }
116   return ret;
117 }
118
119 std::string GetMeshNameWithContext(const std::vector<int>& context)
120 {
121   static const char DFT_MESH_NAME[]="Mesh";
122   if(context.empty())
123     return DFT_MESH_NAME;
124   std::ostringstream oss; oss << DFT_MESH_NAME;
125   for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
126     oss << "_" << *it;
127   return oss.str();
128 }
129
130 DataArrayInt *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
131 {
132   if(!data)
133     throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
134   int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
135   std::size_t nbElts(nbTuples*nbComp);
136   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret(DataArrayInt::New());
137   ret->alloc(nbTuples,nbComp);
138   for(int i=0;i<nbComp;i++)
139     {
140       const char *comp(data->GetComponentName(i));
141       if(comp)
142         ret->setInfoOnComponent(i,comp);
143     }
144   int *ptOut(ret->getPointer());
145   vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
146   if(d0)
147     {
148       const int *pt(d0->GetPointer(0));
149       std::copy(pt,pt+nbElts,ptOut);
150       return ret.retn();
151     }
152   std::ostringstream oss;
153   oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
154   throw MZCException(oss.str());
155 }
156
157 DataArrayDouble *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
158 {
159   if(!data)
160     throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
161   int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
162   std::size_t nbElts(nbTuples*nbComp);
163   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret(DataArrayDouble::New());
164   ret->alloc(nbTuples,nbComp);
165   for(int i=0;i<nbComp;i++)
166     {
167       const char *comp(data->GetComponentName(i));
168       if(comp)
169         ret->setInfoOnComponent(i,comp);
170     }
171   double *ptOut(ret->getPointer());
172   vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
173   if(d0)
174     {
175       const float *pt(d0->GetPointer(0));
176       for(std::size_t i=0;i<nbElts;i++)
177         ptOut[i]=pt[i];
178       return ret.retn();
179     }
180   vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
181   if(d1)
182     {
183       const double *pt(d1->GetPointer(0));
184       std::copy(pt,pt+nbElts,ptOut);
185       return ret.retn();
186     }
187   std::ostringstream oss;
188   oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << typeid(*data).name() << "\" type !";
189   throw MZCException(oss.str());
190 }
191
192 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
193 {
194   if(!data)
195     throw MZCException("ConvertVTKArrayToMCArray : internal error !");
196   vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
197   vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
198   if(d0 || d1)
199     return ConvertVTKArrayToMCArrayDouble(data);
200   vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
201   if(d2)
202     return ConvertVTKArrayToMCArrayInt(data);
203   std::ostringstream oss;
204   oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
205   throw MZCException(oss.str());
206 }
207
208 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
209 {
210   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
211   subMesh->setCoords(coords); subMesh->allocateCells();
212   int nbCells(ca->GetNumberOfCells());
213   if(nbCells==0)
214     return 0;
215   vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
216   const vtkIdType *conn(ca->GetPointer());
217   for(int i=0;i<nbCells;i++)
218     {
219       int sz(*conn++);
220       subMesh->insertNextCell(type,sz,conn);
221       conn+=sz;
222     }
223   return subMesh.retn();
224 }
225
226 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& ids)
227 {
228   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
229   subMesh->setCoords(coords); subMesh->allocateCells();
230   int nbCells(ca->GetNumberOfCells());
231   if(nbCells==0)
232     return 0;
233   vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries());
234   const vtkIdType *conn(ca->GetPointer());
235   ids=DataArrayInt::New() ; ids->alloc(0,1);
236   for(int i=0;i<nbCells;i++)
237     {
238       int sz(*conn++);
239       int nbTri(sz-2);
240       if(nbTri>0)
241         {
242           for(int j=0;j<nbTri;j++,conn++)
243             {
244               subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
245               ids->pushBackSilent(i);
246             }
247         }
248       else
249         {
250           std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
251           throw MZCException(oss.str());
252         }
253       conn+=sz;
254     }
255   return subMesh.retn();
256 }
257
258 class MicroField
259 {
260 public:
261   MicroField(const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh>& m, const std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
262   MicroField(const std::vector< MicroField >& vs);
263   void setNodeFields(const std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> >& nf) { _nodeFs=nf; }
264   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> getMesh() const { return _m; }
265   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > getCellFields() const { return _cellFs; }
266 private:
267   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> _m;
268   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > _cellFs;
269   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > _nodeFs;
270 };
271
272 MicroField::MicroField(const std::vector< MicroField >& vs)
273 {
274   std::size_t sz(vs.size());
275   std::vector<const MEDCouplingUMesh *> vs2(sz);
276   std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<DataArray> > > arrs2(sz);
277   int nbElts(-1);
278   for(std::size_t ii=0;ii<sz;ii++)
279     {
280       vs2[ii]=vs[ii].getMesh();
281       arrs2[ii]=vs[ii].getCellFields();
282       if(nbElts<0)
283         nbElts=arrs2[ii].size();
284       else
285         if(arrs2[ii].size()!=nbElts)
286           throw MZCException("MicroField cstor : internal error !");
287     }
288   _cellFs.resize(nbElts);
289   for(int ii=0;ii<nbElts;ii++)
290     {
291       std::vector<const DataArray *> arrsTmp(sz);
292       for(std::size_t jj=0;jj<sz;jj++)
293         {
294           arrsTmp[jj]=arrs2[jj][ii];
295         }
296       _cellFs[ii]=DataArray::Aggregate(arrsTmp);
297     }
298   _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
299 }
300
301 void AppendMCFieldFrom(ParaMEDMEM::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MEDCouplingAutoRefCountObjectPtr<DataArray> da, const DataArrayInt *n2oPtr)
302 {
303   static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
304   static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
305   if(!da || !mesh || !mfd)
306     throw MZCException("AppendMCFieldFrom : internal error !");
307   MEDFileFields *fs(mfd->getFields());
308   MEDFileMeshes *ms(mfd->getMeshes());
309   if(!fs || !ms)
310     throw MZCException("AppendMCFieldFrom : internal error 2 !");
311   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> dad(ParaMEDMEM::DynamicCast<DataArray,DataArrayDouble>(da));
312   DataArrayDouble *dadPtr(dad);
313   std::string fieldName;
314   if(dadPtr)
315     {
316       fieldName=dadPtr->getName();
317       MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
318       f->setName(fieldName);
319       if(!n2oPtr)
320         f->setArray(dadPtr);
321       else
322         {
323           MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
324           f->setArray(dad2);
325         }
326       f->setMesh(mesh);
327       MEDCouplingAutoRefCountObjectPtr<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
328       MEDCouplingAutoRefCountObjectPtr<MEDFileField1TS> f1ts(MEDFileField1TS::New());
329       f1ts->setFieldNoProfileSBT(f);
330       fmts->pushBackTimeStep(f1ts);
331       fs->pushField(fmts);
332       return ;
333     }
334   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> dai(ParaMEDMEM::DynamicCast<DataArray,DataArrayInt>(da));
335   DataArrayInt *daiPtr(dai);
336   if(daiPtr)
337     {
338       fieldName=daiPtr->getName();
339       if((fieldName!=FAMFIELD_FOR_CELLS || tf!=ParaMEDMEM::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=ParaMEDMEM::ON_NODES))
340         {
341           MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(tf));
342           f->setName(fieldName);
343           f->setMesh(mesh);
344           MEDCouplingAutoRefCountObjectPtr<MEDFileIntFieldMultiTS> fmts(MEDFileIntFieldMultiTS::New());
345           MEDCouplingAutoRefCountObjectPtr<MEDFileIntField1TS> f1ts(MEDFileIntField1TS::New());
346           if(!n2oPtr)
347             f1ts->setFieldNoProfileSBT(f,daiPtr);
348           else
349             {
350               MEDCouplingAutoRefCountObjectPtr<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
351               f1ts->setFieldNoProfileSBT(f,dai2);
352             }
353           fmts->pushBackTimeStep(f1ts);
354           fs->pushField(fmts);
355           return ;
356         }
357       else if(fieldName==FAMFIELD_FOR_CELLS && tf==ParaMEDMEM::ON_CELLS)
358         {
359           MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
360           if(!mm)
361             throw MZCException("AppendMCFieldFrom : internal error 3 !");
362           if(!n2oPtr)
363             mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daiPtr);
364           else
365             {
366               MEDCouplingAutoRefCountObjectPtr<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
367               mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
368             }
369         }
370       else if(fieldName==FAMFIELD_FOR_NODES || tf==ParaMEDMEM::ON_NODES)
371         {
372           MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
373           if(!mm)
374             throw MZCException("AppendMCFieldFrom : internal error 4 !");
375           if(!n2oPtr)
376             mm->setFamilyFieldArr(1,daiPtr);
377           else
378             {
379               MEDCouplingAutoRefCountObjectPtr<DataArrayInt> dai2(daiPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
380               mm->setFamilyFieldArr(1,dai2);
381             }
382         }
383     }
384 }
385
386 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf)
387 {
388   if(!mfd)
389     throw MZCException("PutAtLevelDealOrder : internal error !");
390   MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
391   MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
392   if(!mmu)
393     throw MZCException("PutAtLevelDealOrder : internal error 2 !");
394   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh(mf.getMesh());
395   mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
396   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> o2n(mesh->sortCellsInMEDFileFrmt());
397   const DataArrayInt *o2nPtr(o2n);
398   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> n2o;
399   mmu->setMeshAtLevel(meshDimRel,mesh);
400   const DataArrayInt *n2oPtr(0);
401   if(o2n)
402     {
403       n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
404       n2oPtr=n2o;
405       if(n2oPtr && n2oPtr->isIdentity2(mesh->getNumberOfCells()))
406         n2oPtr=0;
407       if(n2oPtr)
408         mm->setRenumFieldArr(meshDimRel,n2o);
409     }
410   //
411   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cells(mf.getCellFields());
412   for(std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
413     {
414       MEDCouplingAutoRefCountObjectPtr<DataArray> da(*it);
415       AppendMCFieldFrom(ParaMEDMEM::ON_CELLS,mesh,mfd,da,n2oPtr);
416     }
417 }
418
419 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms)
420 {
421   if(!mfd)
422     throw MZCException("AssignSingleGTMeshes : internal error !");
423   MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
424   MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
425   if(!mm)
426     throw MZCException("AssignSingleGTMeshes : internal error 2 !");
427   int meshDim(-std::numeric_limits<int>::max());
428   std::map<int, std::vector< MicroField > > ms2;
429   for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
430     {
431       const MEDCouplingUMesh *elt((*it).getMesh());
432       if(elt)
433         {
434           int myMeshDim(elt->getMeshDimension());
435           meshDim=std::max(meshDim,myMeshDim);
436           ms2[myMeshDim].push_back(*it);
437         }
438     }
439   if(ms2.empty())
440     return ;
441   for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
442     {
443       const std::vector< MicroField >& vs((*it).second);
444       if(vs.size()==1)
445         {
446           PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0]);
447         }
448       else
449         {
450           MicroField merge(vs);
451           PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge);
452         }
453     }
454 }
455
456 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
457 {
458   if(!ds)
459     throw MZCException("BuildCoordsFrom : internal error !");
460   vtkPoints *pts(ds->GetPoints());
461   if(!pts)
462     throw MZCException("BuildCoordsFrom : internal error 2 !");
463   vtkDataArray *data(pts->GetData());
464   if(!data)
465     throw MZCException("BuildCoordsFrom : internal error 3 !");
466   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords(ConvertVTKArrayToMCArrayDouble(data));
467   return coords.retn();
468 }
469
470 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa)
471 {
472   if(!mfd || !dsa)
473     throw MZCException("AddNodeFields : internal error !");
474   MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
475   MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
476   if(!mmu)
477     throw MZCException("AddNodeFields : internal error 2 !");
478   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh(mmu->getMeshAtLevel(0));
479   int nba(dsa->GetNumberOfArrays());
480   for(int i=0;i<nba;i++)
481     {
482       vtkDataArray *arr(dsa->GetArray(i));
483       const char *name(arr->GetName());
484       if(!arr)
485         continue;
486       MEDCouplingAutoRefCountObjectPtr<DataArray> da(ConvertVTKArrayToMCArray(arr));
487       da->setName(name);
488       AppendMCFieldFrom(ParaMEDMEM::ON_NODES,mesh,mfd,da,0);
489     }
490 }
491
492 std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > AddPartFields(const DataArrayInt *part, vtkDataSetAttributes *dsa)
493 {
494   std::vector< MEDCouplingAutoRefCountObjectPtr<DataArray> > ret;
495   if(!dsa)
496     return ret;
497   int nba(dsa->GetNumberOfArrays());
498   for(int i=0;i<nba;i++)
499     {
500       vtkDataArray *arr(dsa->GetArray(i));
501       if(!arr)
502         continue;
503       const char *name(arr->GetName());
504       int nbCompo(arr->GetNumberOfComponents());
505       vtkIdType nbTuples(arr->GetNumberOfTuples());
506       MEDCouplingAutoRefCountObjectPtr<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
507       if(part)
508         mcarr=mcarr->selectByTupleId(part->begin(),part->end());
509       mcarr->setName(name);
510       ret.push_back(mcarr);
511     }
512   return ret;
513 }
514
515 std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
516 {
517   std::vector< MEDCouplingAutoRefCountObjectPtr<DataArray> > ret;
518   if(!dsa)
519     return ret;
520   int nba(dsa->GetNumberOfArrays());
521   for(int i=0;i<nba;i++)
522     {
523       vtkDataArray *arr(dsa->GetArray(i));
524       if(!arr)
525         continue;
526       const char *name(arr->GetName());
527       int nbCompo(arr->GetNumberOfComponents());
528       vtkIdType nbTuples(arr->GetNumberOfTuples());
529       MEDCouplingAutoRefCountObjectPtr<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
530       mcarr=mcarr->selectByTupleId2(bg,end,1);
531       mcarr->setName(name);
532       ret.push_back(mcarr);
533     }
534   return ret;
535 }
536
537 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context)
538 {
539   if(!ds || !ret)
540     throw MZCException("ConvertFromRectilinearGrid : internal error !");
541   //
542   MEDCouplingAutoRefCountObjectPtr<MEDFileMeshes> meshes(MEDFileMeshes::New());
543   ret->setMeshes(meshes);
544   MEDCouplingAutoRefCountObjectPtr<MEDFileFields> fields(MEDFileFields::New());
545   ret->setFields(fields);
546   //
547   MEDCouplingAutoRefCountObjectPtr<MEDFileCMesh> cmesh(MEDFileCMesh::New());
548   meshes->pushMesh(cmesh);
549   MEDCouplingAutoRefCountObjectPtr<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
550   vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
551   if(cx)
552     {
553       MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cx));
554       cmeshmc->setCoordsAt(0,arr);
555     }
556   if(cy)
557     {
558       MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cy));
559       cmeshmc->setCoordsAt(1,arr);
560     }
561   if(cz)
562     {
563       MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDouble(cz));
564       cmeshmc->setCoordsAt(2,arr);
565     }
566   std::string meshName(GetMeshNameWithContext(context));
567   cmeshmc->setName(meshName);
568   cmesh->setMesh(cmeshmc);
569   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
570   for(std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
571     {
572       MEDCouplingAutoRefCountObjectPtr<DataArray> da(*it);
573       AppendMCFieldFrom(ParaMEDMEM::ON_CELLS,cmeshmc,ret,da,0);
574     }
575   std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
576   for(std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
577     {
578       MEDCouplingAutoRefCountObjectPtr<DataArray> da(*it);
579       AppendMCFieldFrom(ParaMEDMEM::ON_NODES,cmeshmc,ret,da,0);
580     }
581 }
582
583 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context)
584 {
585   if(!ds || !ret)
586     throw MZCException("ConvertFromPolyData : internal error !");
587   //
588   MEDCouplingAutoRefCountObjectPtr<MEDFileMeshes> meshes(MEDFileMeshes::New());
589   ret->setMeshes(meshes);
590   MEDCouplingAutoRefCountObjectPtr<MEDFileFields> fields(MEDFileFields::New());
591   ret->setFields(fields);
592   //
593   MEDCouplingAutoRefCountObjectPtr<MEDFileUMesh> umesh(MEDFileUMesh::New());
594   meshes->pushMesh(umesh);
595   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords(BuildCoordsFrom(ds));
596   umesh->setCoords(coords);
597   umesh->setName(GetMeshNameWithContext(context));
598   //
599   int offset(0);
600   std::vector< MicroField > ms;
601   vtkCellArray *cd(ds->GetVerts());
602   if(cd)
603     {
604       MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
605       if((const MEDCouplingUMesh *)subMesh)
606         {
607           std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
608           offset+=subMesh->getNumberOfCells();
609           ms.push_back(MicroField(subMesh,cellFs));
610         }
611     }
612   vtkCellArray *cc(ds->GetLines());
613   if(cc)
614     {
615       MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2));
616       if((const MEDCouplingUMesh *)subMesh)
617         {
618           std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
619           offset+=subMesh->getNumberOfCells();
620           ms.push_back(MicroField(subMesh,cellFs));
621         }
622     }
623   vtkCellArray *cb(ds->GetPolys());
624   if(cb)
625     {
626       MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
627       if((const MEDCouplingUMesh *)subMesh)
628         {
629           std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
630           offset+=subMesh->getNumberOfCells();
631           ms.push_back(MicroField(subMesh,cellFs));
632         }
633     }
634   vtkCellArray *ca(ds->GetStrips());
635   if(ca)
636     {
637       MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids;
638       MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
639       if((const MEDCouplingUMesh *)subMesh)
640         {
641           std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
642           offset+=subMesh->getNumberOfCells();
643           ms.push_back(MicroField(subMesh,cellFs));
644         }
645     }
646   AssignSingleGTMeshes(ret,ms);
647   AddNodeFields(ret,ds->GetPointData());
648 }
649
650 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context)
651 {
652   if(!ds || !ret)
653     throw MZCException("ConvertFromUnstructuredGrid : internal error !");
654   //
655   MEDCouplingAutoRefCountObjectPtr<MEDFileMeshes> meshes(MEDFileMeshes::New());
656   ret->setMeshes(meshes);
657   MEDCouplingAutoRefCountObjectPtr<MEDFileFields> fields(MEDFileFields::New());
658   ret->setFields(fields);
659   //
660   MEDCouplingAutoRefCountObjectPtr<MEDFileUMesh> umesh(MEDFileUMesh::New());
661   meshes->pushMesh(umesh);
662   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords(BuildCoordsFrom(ds));
663   umesh->setCoords(coords);
664   umesh->setName(GetMeshNameWithContext(context));
665   vtkIdType nbCells(ds->GetNumberOfCells());
666   vtkCellArray *ca(ds->GetCells());
667   if(!ca)
668     return ;
669   vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries());
670   vtkIdType *caPtr(ca->GetPointer());
671   vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
672   if(!ct)
673     throw MZCException("ConvertFromUnstructuredGrid : internal error");
674   vtkIdTypeArray *cla(ds->GetCellLocationsArray());
675   const vtkIdType *claPtr(cla->GetPointer(0));
676   if(!cla)
677     throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
678   const unsigned char *ctPtr(ct->GetPointer(0));
679   std::map<int,int> m(ComputeMapOfType());
680   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> lev(DataArrayInt::New()) ;  lev->alloc(nbCells,1);
681   int *levPtr(lev->getPointer());
682   for(vtkIdType i=0;i<nbCells;i++)
683     {
684       std::map<int,int>::iterator it(m.find(ctPtr[i]));
685       if(it!=m.end())
686         {
687           const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
688           levPtr[i]=cm.getDimension();
689         }
690       else
691         {
692           std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
693           throw MZCException(oss.str());
694         }
695     }
696   int dummy(0);
697   int meshDim(lev->getMaxValue(dummy));
698   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> levs(lev->getDifferentValues());
699   std::vector< MicroField > ms;
700   vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations());
701   for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
702     {
703       MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
704       m0->setCoords(coords); m0->allocateCells();
705       MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cellIdsCurLev(lev->getIdsEqual(*curLev));
706       for(const int *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
707         {
708           std::map<int,int>::iterator it(m.find(ctPtr[*cellId]));
709           vtkIdType offset(claPtr[*cellId]);
710           vtkIdType sz(caPtr[offset]);
711           INTERP_KERNEL::NormalizedCellType ct((INTERP_KERNEL::NormalizedCellType)(*it).second);
712           if(ct!=INTERP_KERNEL::NORM_POLYHED)
713             m0->insertNextCell(ct,sz,caPtr+offset+1);
714           else
715             {
716               if(!faces || !faceLoc)
717                 throw MZCException("ConvertFromUnstructuredGrid : faces are expected when there are polyhedra !");
718               const vtkIdType *facPtr(faces->GetPointer(0)),*facLocPtr(faceLoc->GetPointer(0));
719               std::vector<int> conn;
720               int off0(facLocPtr[*cellId]);
721               int nbOfFaces(facPtr[off0++]);
722               for(int k=0;k<nbOfFaces;k++)
723                 {
724                   int nbOfNodesInFace(facPtr[off0++]);
725                   std::copy(facPtr+off0,facPtr+off0+nbOfNodesInFace,std::back_inserter(conn));
726                   off0+=nbOfNodesInFace;
727                   if(k<nbOfFaces-1)
728                     conn.push_back(-1);
729                 }
730               m0->insertNextCell(ct,conn.size(),&conn[0]);
731             }
732         }
733       std::vector<MEDCouplingAutoRefCountObjectPtr<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
734       ms.push_back(MicroField(m0,cellFs));
735     }
736   AssignSingleGTMeshes(ret,ms);
737   AddNodeFields(ret,ds->GetPointData());
738 }
739
740 ///////////////////
741
742 vtkMEDWriter::vtkMEDWriter():WriteAllTimeSteps(0),FileName(0),IsTouched(false)
743 {
744 }
745
746 vtkMEDWriter::~vtkMEDWriter()
747 {
748 }
749
750 vtkInformationDataObjectMetaDataKey *GetMEDReaderMetaDataIfAny()
751 {
752   static const char ZE_KEY[]="vtkMEDReader::META_DATA";
753   ParaMEDMEM::GlobalDict *gd(ParaMEDMEM::GlobalDict::GetInstance());
754   if(!gd->hasKey(ZE_KEY))
755     return 0;
756   std::string ptSt(gd->value(ZE_KEY));
757   void *pt(0);
758   std::istringstream iss(ptSt); iss >> pt;
759   return reinterpret_cast<vtkInformationDataObjectMetaDataKey *>(pt);
760 }
761
762 bool IsInformationOK(vtkInformation *info)
763 {
764   vtkInformationDataObjectMetaDataKey *key(GetMEDReaderMetaDataIfAny());
765   if(!key)
766     return false;
767   // Check the information contain meta data key
768   if(!info->Has(key))
769     return false;
770   // Recover Meta Data
771   vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
772   if(!sil)
773     return false;
774   int idNames(0);
775   vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
776   vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
777   if(!verticesNames2)
778     return false;
779   for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
780     {
781       vtkStdString &st(verticesNames2->GetValue(i));
782       if(st=="MeshesFamsGrps")
783         return true;
784     }
785   return false;
786 }
787
788 class Grp
789 {
790 public:
791   Grp(const std::string& name):_name(name) { }
792   void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
793   std::string getName() const { return _name; }
794   std::vector<std::string> getFamilies() const { return _fams; }
795 private:
796   std::string _name;
797   std::vector<std::string> _fams;
798 };
799
800 class Fam
801 {
802 public:
803   Fam(const std::string& name)
804   {
805     static const char ZE_SEP[]="@@][@@";
806     std::size_t pos(name.find(ZE_SEP));
807     std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
808     std::istringstream iss(name1);
809     iss >> _id;
810     _name=name0;
811   }
812   std::string getName() const { return _name; }
813   int getID() const { return _id; }
814 private:
815   std::string _name;
816   int _id;
817 };
818
819 void LoadFamGrpMapInfo(vtkMutableDirectedGraph *sil, std::string& meshName, std::vector<Grp>& groups, std::vector<Fam>& fams)
820 {
821   if(!sil)
822     throw MZCException("LoadFamGrpMapInfo : internal error !");
823   int idNames(0);
824   vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
825   vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
826   vtkIdType id0;
827   bool found(false);
828   for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
829     {
830       vtkStdString &st(verticesNames2->GetValue(i));
831       if(st=="MeshesFamsGrps")
832         {
833           id0=i;
834           found=true;
835         }
836     }
837   if(!found)
838     throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
839   vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
840   sil->GetAdjacentVertices(id0,it0);
841   int kk(0),ll(0);
842   while(it0->HasNext())
843     {
844       vtkIdType id1(it0->Next());
845       std::string mName((const char *)verticesNames2->GetValue(id1));
846       meshName=mName;
847       vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
848       sil->GetAdjacentVertices(id1,it1);
849       vtkIdType idZeGrps(it1->Next());//zeGroups
850       vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
851       sil->GetAdjacentVertices(idZeGrps,itGrps);
852       while(itGrps->HasNext())
853         {
854           vtkIdType idg(itGrps->Next());
855           Grp grp((const char *)verticesNames2->GetValue(idg));
856           vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
857           sil->GetAdjacentVertices(idg,itGrps2);
858           std::vector<std::string> famsOnGroup;
859           while(itGrps2->HasNext())
860             {
861               vtkIdType idgf(itGrps2->Next());
862               famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
863             }
864           grp.setFamilies(famsOnGroup);
865           itGrps2->Delete();
866           groups.push_back(grp);
867         }
868       itGrps->Delete();
869       vtkIdType idZeFams(it1->Next());//zeFams
870       it1->Delete();
871       vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
872       sil->GetAdjacentVertices(idZeFams,itFams);
873       while(itFams->HasNext())
874         {
875           vtkIdType idf(itFams->Next());
876           Fam fam((const char *)verticesNames2->GetValue(idf));
877           fams.push_back(fam);
878         }
879       itFams->Delete();
880     }
881   it0->Delete();
882 }
883
884 int vtkMEDWriter::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
885
886   //std::cerr << "########################################## vtkMEDWriter::RequestInformation ########################################## " << (const char *) this->FileName << std::endl;
887   return 1;
888 }
889
890 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context)
891 {
892   if(!ds || !mfd)
893     throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
894   vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
895   vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
896   vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
897   if(ds2)
898     {
899       ConvertFromPolyData(mfd,ds2,context);
900     }
901   else if(ds3)
902     {
903       ConvertFromUnstructuredGrid(mfd,ds3,context);
904     }
905   else if(ds4)
906     {
907       ConvertFromRectilinearGrid(mfd,ds4,context);
908     }
909   else
910     throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
911 }
912
913 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context)
914 {
915   if(!ds || !mfd)
916     throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
917   int nbBlocks(ds->GetNumberOfBlocks());
918   if(nbBlocks==1 && context.empty())
919     {
920       vtkDataObject *uniqueElt(ds->GetBlock(0));
921       if(!uniqueElt)
922         throw MZCException("Unique elt in multiblock is NULL !");
923       vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
924       if(uniqueEltc)
925         {
926           WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context);
927           return ;
928         }
929     }
930   for(int i=0;i<nbBlocks;i++)
931     {
932       vtkDataObject *elt(ds->GetBlock(i));
933       std::vector<int> context2;
934       context2.push_back(i);
935       if(!elt)
936         {
937           std::ostringstream oss; oss << "In context ";
938           std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
939           oss << " at pos #" << i << " elt is NULL !";
940           throw MZCException(oss.str());
941         }
942       vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
943       if(elt1)
944         {
945           WriteMEDFileFromVTKDataSet(mfd,elt1,context);
946           continue;
947         }
948       vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
949       if(elt2)
950         {
951           WriteMEDFileFromVTKMultiBlock(mfd,elt2,context);
952           continue;
953         }
954       std::ostringstream oss; oss << "In context ";
955       std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
956       oss << " at pos #" << i << " elt not recognized data type !";
957       throw MZCException(oss.str());
958     }
959 }
960
961 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input)
962 {
963   if(!input || !mfd)
964     throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
965   std::vector<int> context;
966   vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
967   if(input1)
968     {
969       WriteMEDFileFromVTKDataSet(mfd,input1,context);
970       return ;
971     }
972   vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
973   if(input2)
974     {
975       WriteMEDFileFromVTKMultiBlock(mfd,input2,context);
976       return ;
977     }
978   throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
979 }
980
981 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
982 {
983   if(!mfd)
984     return ;
985   if(meshName.empty())
986     return ;
987   MEDFileMeshes *meshes(mfd->getMeshes());
988   if(!meshes)
989     return ;
990   if(meshes->getNumberOfMeshes()!=1)
991     return ;
992   MEDFileMesh *mm(meshes->getMeshAtPos(0));
993   if(!mm)
994     return ;
995   mm->setName(meshName);
996   for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
997     mm->setFamilyId((*it).getName(),(*it).getID());
998   for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
999     mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
1000   MEDFileFields *fields(mfd->getFields());
1001   if(!fields)
1002     return ;
1003   for(int i=0;i<fields->getNumberOfFields();i++)
1004     {
1005       MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
1006       if(!fmts)
1007         continue;
1008       fmts->setMeshName(meshName);
1009     }
1010 }
1011
1012 int vtkMEDWriter::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
1013 {
1014   //std::cerr << "########################################## vtkMEDWriter::RequestData        ########################################## " << (const char *) this->FileName << std::endl;
1015   try
1016     {
1017       vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
1018       std::string meshName;
1019       std::vector<Grp> groups; std::vector<Fam> fams;
1020       if(IsInformationOK(inputInfo))
1021         {
1022           vtkMutableDirectedGraph *famGrpGraph(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
1023           LoadFamGrpMapInfo(famGrpGraph,meshName,groups,fams);
1024         }
1025       vtkInformation *outInfo(outputVector->GetInformationObject(0));
1026       vtkDataObject *input(vtkDataObject::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
1027       if(!input)
1028         throw MZCException("Not recognized data object in input of the MEDWriter ! Maybe not implemented yet !");
1029       MEDCouplingAutoRefCountObjectPtr<MEDFileData> mfd(MEDFileData::New());
1030       WriteMEDFileFromVTKGDS(mfd,input);
1031       PutFamGrpInfoIfAny(mfd,meshName,groups,fams);
1032       mfd->write(this->FileName,this->IsTouched?0:2); this->IsTouched=true;
1033       outInfo->Set(vtkDataObject::DATA_OBJECT(),input);
1034     }
1035   catch(MZCException& e)
1036     {
1037       std::ostringstream oss;
1038       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;
1039       if(this->HasObserver("ErrorEvent") )
1040         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
1041       else
1042         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
1043       vtkObject::BreakOnError();
1044       return 0;
1045     }
1046   return 1;
1047 }
1048
1049 void vtkMEDWriter::PrintSelf(ostream& os, vtkIndent indent)
1050 {
1051   this->Superclass::PrintSelf(os, indent);
1052 }
1053
1054 int vtkMEDWriter::Write()
1055 {
1056   this->Update();
1057   return 0;
1058 }