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