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