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