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