Salome HOME
Try to solve MEDWRITER_0 test on PV590
[tools/paravisaddons_common.git] / src / RateOfFlowThroughSection / plugin / RateOfFlowThroughSectionModule / VTKToMEDMem.cxx
1 // Copyright (C) 2017-2021  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 "VTKToMEDMem.h"
22
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkIntArray.h"
25 #include "vtkLongArray.h"
26 #include "vtkLongLongArray.h"
27 #include "vtkCellData.h"
28 #include "vtkPointData.h"
29 #include "vtkFloatArray.h"
30 #include "vtkCellArray.h"
31
32 #include "vtkStreamingDemandDrivenPipeline.h"
33 #include "vtkInformationDataObjectMetaDataKey.h"
34 #include "vtkUnstructuredGrid.h"
35 #include "vtkMultiBlockDataSet.h"
36 #include "vtkRectilinearGrid.h"
37 #include "vtkInformationStringKey.h"
38 #include "vtkAlgorithmOutput.h"
39 #include "vtkObjectFactory.h"
40 #include "vtkMutableDirectedGraph.h"
41 #include "vtkMultiBlockDataSet.h"
42 #include "vtkPolyData.h"
43 #include "vtkDataSet.h"
44 #include "vtkInformationVector.h"
45 #include "vtkInformation.h"
46 #include "vtkDataArraySelection.h"
47 #include "vtkTimeStamp.h"
48 #include "vtkInEdgeIterator.h"
49 #include "vtkInformationDataObjectKey.h"
50 #include "vtkExecutive.h"
51 #include "vtkVariantArray.h"
52 #include "vtkStringArray.h"
53 #include "vtkDoubleArray.h"
54 #include "vtkCharArray.h"
55 #include "vtkUnsignedCharArray.h"
56 #include "vtkDataSetAttributes.h"
57 #include "vtkDemandDrivenPipeline.h"
58 #include "vtkDataObjectTreeIterator.h"
59 #include "vtkWarpScalar.h"
60
61 #include <map>
62 #include <deque>
63 #include <sstream>
64 #include <cstring>
65
66 using VTKToMEDMem::Grp;
67 using VTKToMEDMem::Fam;
68
69 using MEDCoupling::MEDFileData;
70 using MEDCoupling::MEDFileMesh;
71 using MEDCoupling::MEDFileCMesh;
72 using MEDCoupling::MEDFileUMesh;
73 using MEDCoupling::MEDFileFields;
74 using MEDCoupling::MEDFileMeshes;
75
76 using MEDCoupling::MEDFileInt32Field1TS;
77 using MEDCoupling::MEDFileInt64Field1TS;
78 using MEDCoupling::MEDFileField1TS;
79 using MEDCoupling::MEDFileIntFieldMultiTS;
80 using MEDCoupling::MEDFileFieldMultiTS;
81 using MEDCoupling::MEDFileAnyTypeFieldMultiTS;
82 using MEDCoupling::DataArray;
83 using MEDCoupling::DataArrayInt32;
84 using MEDCoupling::DataArrayInt64;
85 using MEDCoupling::DataArrayFloat;
86 using MEDCoupling::DataArrayDouble;
87 using MEDCoupling::MEDCouplingMesh;
88 using MEDCoupling::MEDCouplingUMesh;
89 using MEDCoupling::MEDCouplingCMesh;
90 using MEDCoupling::MEDCouplingFieldDouble;
91 using MEDCoupling::MEDCouplingFieldFloat;
92 using MEDCoupling::MEDCouplingFieldInt;
93 using MEDCoupling::MEDCouplingFieldInt64;
94 using MEDCoupling::MCAuto;
95 using MEDCoupling::Traits;
96 using MEDCoupling::MLFieldTraits;
97
98 ///////////////////
99
100 Fam::Fam(const std::string& name)
101 {
102   static const char ZE_SEP[]="@@][@@";
103   std::size_t pos(name.find(ZE_SEP));
104   std::string name0(name.substr(0,pos)),name1(name.substr(pos+strlen(ZE_SEP)));
105   std::istringstream iss(name1);
106   iss >> _id;
107   _name=name0;
108 }
109
110 ///////////////////
111
112 #include "VTKMEDTraits.hxx"
113
114 std::map<int,int> ComputeMapOfType()
115 {
116   std::map<int,int> ret;
117   int nbOfTypesInMC(sizeof(MEDCOUPLING2VTKTYPETRADUCER)/sizeof( decltype(MEDCOUPLING2VTKTYPETRADUCER[0]) ));
118   for(int i=0;i<nbOfTypesInMC;i++)
119     {
120       auto vtkId(MEDCOUPLING2VTKTYPETRADUCER[i]);
121       if(vtkId!=MEDCOUPLING2VTKTYPETRADUCER_NONE)
122         ret[vtkId]=i;
123     }
124   return ret;
125 }
126
127 std::string GetMeshNameWithContext(const std::vector<int>& context)
128 {
129   static const char DFT_MESH_NAME[]="Mesh";
130   if(context.empty())
131     return DFT_MESH_NAME;
132   std::ostringstream oss; oss << DFT_MESH_NAME;
133   for(std::vector<int>::const_iterator it=context.begin();it!=context.end();it++)
134     oss << "_" << *it;
135   return oss.str();
136 }
137
138 DataArrayIdType *ConvertVTKArrayToMCArrayInt(vtkDataArray *data)
139 {
140   if(!data)
141     throw MZCException("ConvertVTKArrayToMCArrayInt : internal error !");
142   int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
143   std::size_t nbElts(nbTuples*nbComp);
144   MCAuto<DataArrayIdType> ret(DataArrayIdType::New());
145   ret->alloc(nbTuples,nbComp);
146   for(int i=0;i<nbComp;i++)
147     {
148       const char *comp(data->GetComponentName(i));
149       if(comp)
150         ret->setInfoOnComponent(i,comp);
151     }
152   mcIdType *ptOut(ret->getPointer());
153   vtkIntArray *d0(vtkIntArray::SafeDownCast(data));
154   if(d0)
155     {
156       const int *pt(d0->GetPointer(0));
157       std::copy(pt,pt+nbElts,ptOut);
158       return ret.retn();
159     }
160   vtkLongArray *d1(vtkLongArray::SafeDownCast(data));
161   if(d1)
162     {
163       const long *pt(d1->GetPointer(0));
164       std::copy(pt,pt+nbElts,ptOut);
165       return ret.retn();
166     }
167   vtkLongLongArray *d1l(vtkLongLongArray::SafeDownCast(data));
168   if(d1l)
169     {
170       const long long *pt(d1l->GetPointer(0));
171       std::copy(pt,pt+nbElts,ptOut);
172       return ret.retn();
173     }
174   vtkIdTypeArray *d3(vtkIdTypeArray::SafeDownCast(data));
175   if(d3)
176     {
177       const vtkIdType *pt(d3->GetPointer(0));
178       std::copy(pt,pt+nbElts,ptOut);
179       return ret.retn();
180     }
181   vtkUnsignedCharArray *d2(vtkUnsignedCharArray::SafeDownCast(data));
182   if(d2)
183     {
184       const unsigned char *pt(d2->GetPointer(0));
185       std::copy(pt,pt+nbElts,ptOut);
186       return ret.retn();
187     }
188   std::ostringstream oss;
189   oss << "ConvertVTKArrayToMCArrayInt : unrecognized array \"" << typeid(*data).name() << "\" type !";
190   throw MZCException(oss.str());
191 }
192
193 template<class T>
194 typename Traits<T>::ArrayType *ConvertVTKArrayToMCArrayDouble(vtkDataArray *data)
195 {
196   if(!data)
197     throw MZCException("ConvertVTKArrayToMCArrayDouble : internal error !");
198   int nbTuples(data->GetNumberOfTuples()),nbComp(data->GetNumberOfComponents());
199   std::size_t nbElts(nbTuples*nbComp);
200   MCAuto< typename Traits<T>::ArrayType > ret(Traits<T>::ArrayType::New());
201   ret->alloc(nbTuples,nbComp);
202   for(int i=0;i<nbComp;i++)
203     {
204       const char *comp(data->GetComponentName(i));
205       if(comp)
206         ret->setInfoOnComponent(i,comp);
207       else
208         {
209           if(nbComp>1 && nbComp<=3)
210             {
211               char tmp[2];
212               tmp[0]=(char)('X'+i); tmp[1]='\0';
213               ret->setInfoOnComponent(i,tmp);
214             }
215         }
216     }
217   T *ptOut(ret->getPointer());
218   typename MEDFileVTKTraits<T>::VtkType *d0(MEDFileVTKTraits<T>::VtkType::SafeDownCast(data));
219   if(d0)
220     {
221       const T *pt(d0->GetPointer(0));
222       for(std::size_t i=0;i<nbElts;i++)
223         ptOut[i]=(T)pt[i];
224       return ret.retn();
225     }
226   std::ostringstream oss;
227   oss << "ConvertVTKArrayToMCArrayDouble : unrecognized array \"" << data->GetClassName() << "\" type !";
228   throw MZCException(oss.str());
229 }
230
231 DataArrayDouble *ConvertVTKArrayToMCArrayDoubleForced(vtkDataArray *data)
232 {
233   if(!data)
234     throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : internal error 0 !");
235   vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
236   if(d0)
237     {
238       MCAuto<DataArrayFloat> ret(ConvertVTKArrayToMCArrayDouble<float>(data));
239       MCAuto<DataArrayDouble> ret2(ret->convertToDblArr());
240       return ret2.retn();
241     }
242   vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
243   if(d1)
244     return ConvertVTKArrayToMCArrayDouble<double>(data);
245   throw MZCException("ConvertVTKArrayToMCArrayDoubleForced : unrecognized type of data for double !");
246 }
247
248 DataArray *ConvertVTKArrayToMCArray(vtkDataArray *data)
249 {
250   if(!data)
251     throw MZCException("ConvertVTKArrayToMCArray : internal error !");
252   vtkFloatArray *d0(vtkFloatArray::SafeDownCast(data));
253   if(d0)
254     return ConvertVTKArrayToMCArrayDouble<float>(data);
255   vtkDoubleArray *d1(vtkDoubleArray::SafeDownCast(data));
256   if(d1)
257     return ConvertVTKArrayToMCArrayDouble<double>(data);
258   vtkIntArray *d2(vtkIntArray::SafeDownCast(data));
259   vtkLongArray *d3(vtkLongArray::SafeDownCast(data));
260   vtkUnsignedCharArray *d4(vtkUnsignedCharArray::SafeDownCast(data));
261   vtkIdTypeArray *d5(vtkIdTypeArray::SafeDownCast(data));
262   if(d2 || d3 || d4 || d5)
263     return ConvertVTKArrayToMCArrayInt(data);
264   std::ostringstream oss;
265   oss << "ConvertVTKArrayToMCArray : unrecognized array \"" << typeid(*data).name() << "\" type !";
266   throw MZCException(oss.str());
267 }
268
269 MEDCouplingUMesh *BuildMeshFromCellArray(vtkCellArray *ca, DataArrayDouble *coords, int meshDim, INTERP_KERNEL::NormalizedCellType type)
270 {
271   MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",meshDim));
272   subMesh->setCoords(coords); subMesh->allocateCells();
273   int nbCells(ca->GetNumberOfCells());
274   if(nbCells==0)
275     return 0;
276   //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
277   const vtkIdType *conn(ca->GetData()->GetPointer(0));
278   for(int i=0;i<nbCells;i++)
279     {
280       mcIdType sz(ToIdType(*conn++));
281       std::vector<mcIdType> conn2(sz);
282       for(int jj=0;jj<sz;jj++)
283         conn2[jj]=ToIdType(conn[jj]);
284       subMesh->insertNextCell(type,sz,&conn2[0]);
285       conn+=sz;
286     }
287   return subMesh.retn();
288 }
289
290 MEDCouplingUMesh *BuildMeshFromCellArrayTriangleStrip(vtkCellArray *ca, DataArrayDouble *coords, MCAuto<DataArrayIdType>& ids)
291 {
292   MCAuto<MEDCouplingUMesh> subMesh(MEDCouplingUMesh::New("",2));
293   subMesh->setCoords(coords); subMesh->allocateCells();
294   int nbCells(ca->GetNumberOfCells());
295   if(nbCells==0)
296     return 0;
297   //vtkIdType nbEntries(ca->GetNumberOfConnectivityEntries()); // todo: unused
298   const vtkIdType *conn(ca->GetData()->GetPointer(0));
299   ids=DataArrayIdType::New() ; ids->alloc(0,1);
300   for(int i=0;i<nbCells;i++)
301     {
302       int sz(*conn++);
303       int nbTri(sz-2);
304       if(nbTri>0)
305         {
306           for(int j=0;j<nbTri;j++,conn++)
307             {
308               mcIdType conn2[3]; conn2[0]=conn[0] ; conn2[1]=conn[1] ; conn2[2]=conn[2];
309               subMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn2);
310               ids->pushBackSilent(i);
311             }
312         }
313       else
314         {
315           std::ostringstream oss; oss << "BuildMeshFromCellArrayTriangleStrip : on cell #" << i << " the triangle stip looks bab !";
316           throw MZCException(oss.str());
317         }
318       conn+=sz;
319     }
320   return subMesh.retn();
321 }
322
323 class MicroField
324 {
325 public:
326   MicroField(const MCAuto<MEDCouplingUMesh>& m, const std::vector<MCAuto<DataArray> >& cellFs):_m(m),_cellFs(cellFs) { }
327   MicroField(const std::vector< MicroField >& vs);
328   void setNodeFields(const std::vector<MCAuto<DataArray> >& nf) { _nodeFs=nf; }
329   MCAuto<MEDCouplingUMesh> getMesh() const { return _m; }
330   std::vector<MCAuto<DataArray> > getCellFields() const { return _cellFs; }
331 private:
332   MCAuto<MEDCouplingUMesh> _m;
333   std::vector<MCAuto<DataArray> > _cellFs;
334   std::vector<MCAuto<DataArray> > _nodeFs;
335 };
336
337 MicroField::MicroField(const std::vector< MicroField >& vs)
338 {
339   std::size_t sz(vs.size());
340   std::vector<const MEDCouplingUMesh *> vs2(sz);
341   std::vector< std::vector< MCAuto<DataArray> > > arrs2(sz);
342   int nbElts(-1);
343   for(std::size_t ii=0;ii<sz;ii++)
344     {
345       vs2[ii]=vs[ii].getMesh();
346       arrs2[ii]=vs[ii].getCellFields();
347       if(nbElts<0)
348         nbElts=(int)arrs2[ii].size();
349       else
350         if((int)arrs2[ii].size()!=nbElts)
351           throw MZCException("MicroField cstor : internal error !");
352     }
353   _cellFs.resize(nbElts);
354   for(int ii=0;ii<nbElts;ii++)
355     {
356       std::vector<const DataArray *> arrsTmp(sz);
357       for(std::size_t jj=0;jj<sz;jj++)
358         {
359           arrsTmp[jj]=arrs2[jj][ii];
360         }
361       _cellFs[ii]=DataArray::Aggregate(arrsTmp);
362     }
363   _m=MEDCouplingUMesh::MergeUMeshesOnSameCoords(vs2);
364 }
365
366 template<class T>
367 void AppendToFields(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, const DataArrayIdType *n2oPtr, typename MEDFileVTKTraits<T>::MCType *dadPtr, MEDFileFields *fs, double timeStep, int tsId)
368 {
369   std::string fieldName(dadPtr->getName());
370   MCAuto< typename Traits<T>::FieldType > f(Traits<T>::FieldType::New(tf));
371   f->setTime(timeStep,tsId,0);
372   {
373     std::string fieldNameForChuckNorris(MEDCoupling::MEDFileAnyTypeField1TSWithoutSDA::FieldNameToMEDFileConvention(fieldName));
374     f->setName(fieldNameForChuckNorris);
375   }
376   if(!n2oPtr)
377     f->setArray(dadPtr);
378   else
379     {
380       MCAuto< typename Traits<T>::ArrayType > dad2(dadPtr->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
381       f->setArray(dad2);
382     }
383   f->setMesh(mesh);
384   MCAuto< typename MLFieldTraits<T>::FMTSType > fmts(MLFieldTraits<T>::FMTSType::New());
385   MCAuto< typename MLFieldTraits<T>::F1TSType > f1ts(MLFieldTraits<T>::F1TSType::New());
386   f1ts->setFieldNoProfileSBT(f);
387   fmts->pushBackTimeStep(f1ts);
388   fs->pushField(fmts);
389 }
390
391 void AppendMCFieldFrom(MEDCoupling::TypeOfField tf, MEDCouplingMesh *mesh, MEDFileData *mfd, MCAuto<DataArray> da, const DataArrayIdType *n2oPtr, double timeStep, int tsId)
392 {
393   static const char FAMFIELD_FOR_CELLS[]="FamilyIdCell";
394   static const char FAMFIELD_FOR_NODES[]="FamilyIdNode";
395   if(!da || !mesh || !mfd)
396     throw MZCException("AppendMCFieldFrom : internal error !");
397   MEDFileFields *fs(mfd->getFields());
398   MEDFileMeshes *ms(mfd->getMeshes());
399   if(!fs || !ms)
400     throw MZCException("AppendMCFieldFrom : internal error 2 !");
401   MCAuto<DataArrayDouble> dad(MEDCoupling::DynamicCast<DataArray,DataArrayDouble>(da));
402   if(dad.isNotNull())
403     {
404       AppendToFields<double>(tf,mesh,n2oPtr,dad,fs,timeStep,tsId);
405       return ;
406     }
407   MCAuto<DataArrayFloat> daf(MEDCoupling::DynamicCast<DataArray,DataArrayFloat>(da));
408   if(daf.isNotNull())
409     {
410       AppendToFields<float>(tf,mesh,n2oPtr,daf,fs,timeStep,tsId);
411       return ;
412     }
413   MCAuto<DataArrayInt> dai(MEDCoupling::DynamicCast<DataArray,DataArrayInt>(da));
414   MCAuto<DataArrayIdType> daId(MEDCoupling::DynamicCast<DataArray,DataArrayIdType>(da));
415   if(dai.isNotNull() || daId.isNotNull())
416     {
417       std::string fieldName(da->getName());
418       if((fieldName!=FAMFIELD_FOR_CELLS || tf!=MEDCoupling::ON_CELLS) && (fieldName!=FAMFIELD_FOR_NODES || tf!=MEDCoupling::ON_NODES))
419         {
420           if(!dai)
421             AppendToFields<mcIdType>(tf,mesh,n2oPtr,daId,fs,timeStep,tsId);
422           else
423             AppendToFields<int>(tf,mesh,n2oPtr,dai,fs,timeStep,tsId);
424           return ;
425         }
426       else if(fieldName==FAMFIELD_FOR_CELLS && tf==MEDCoupling::ON_CELLS)
427         {
428           MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
429           if(!mm)
430             throw MZCException("AppendMCFieldFrom : internal error 3 !");
431           if(!daId)
432             throw MZCException("AppendMCFieldFrom : internal error 3 (not mcIdType) !");
433           if(!n2oPtr)
434             mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),daId);
435           else
436             {
437               MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
438               mm->setFamilyFieldArr(mesh->getMeshDimension()-mm->getMeshDimension(),dai2);
439             }
440         }
441       else if(fieldName==FAMFIELD_FOR_NODES || tf==MEDCoupling::ON_NODES)
442         {
443           MEDFileMesh *mm(ms->getMeshWithName(mesh->getName()));
444           if(!mm)
445             throw MZCException("AppendMCFieldFrom : internal error 4 !");
446           if(!daId)
447             throw MZCException("AppendMCFieldFrom : internal error 4 (not mcIdType) !");
448           if(!n2oPtr)
449             mm->setFamilyFieldArr(1,daId);
450           else
451             {
452               MCAuto<DataArrayIdType> dai2(daId->selectByTupleId(n2oPtr->begin(),n2oPtr->end()));
453               mm->setFamilyFieldArr(1,dai2);
454             }
455         }
456     }
457 }
458
459 void PutAtLevelDealOrder(MEDFileData *mfd, int meshDimRel, const MicroField& mf, double timeStep, int tsId)
460 {
461   if(!mfd)
462     throw MZCException("PutAtLevelDealOrder : internal error !");
463   MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
464   MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
465   if(!mmu)
466     throw MZCException("PutAtLevelDealOrder : internal error 2 !");
467   MCAuto<MEDCouplingUMesh> mesh(mf.getMesh());
468   mesh->setName(mfd->getMeshes()->getMeshAtPos(0)->getName());
469   MCAuto<DataArrayIdType> o2n(mesh->sortCellsInMEDFileFrmt());
470   //const DataArrayIdType *o2nPtr(o2n); // todo: unused
471   MCAuto<DataArrayIdType> n2o;
472   mmu->setMeshAtLevel(meshDimRel,mesh);
473   const DataArrayIdType *n2oPtr(0);
474   if(o2n)
475     {
476       n2o=o2n->invertArrayO2N2N2O(mesh->getNumberOfCells());
477       n2oPtr=n2o;
478       if(n2oPtr && n2oPtr->isIota(mesh->getNumberOfCells()))
479         n2oPtr=0;
480       if(n2oPtr)
481         mm->setRenumFieldArr(meshDimRel,n2o);
482     }
483   //
484   std::vector<MCAuto<DataArray> > cells(mf.getCellFields());
485   for(std::vector<MCAuto<DataArray> >::const_iterator it=cells.begin();it!=cells.end();it++)
486     {
487       MCAuto<DataArray> da(*it);
488       AppendMCFieldFrom(MEDCoupling::ON_CELLS,mesh,mfd,da,n2oPtr,timeStep,tsId);
489     }
490 }
491
492 void AssignSingleGTMeshes(MEDFileData *mfd, const std::vector< MicroField >& ms, double timeStep, int tsId)
493 {
494   if(!mfd)
495     throw MZCException("AssignSingleGTMeshes : internal error !");
496   MEDFileMesh *mm0(mfd->getMeshes()->getMeshAtPos(0));
497   MEDFileUMesh *mm(dynamic_cast<MEDFileUMesh *>(mm0));
498   if(!mm)
499     throw MZCException("AssignSingleGTMeshes : internal error 2 !");
500   int meshDim(-std::numeric_limits<int>::max());
501   std::map<int, std::vector< MicroField > > ms2;
502   for(std::vector< MicroField >::const_iterator it=ms.begin();it!=ms.end();it++)
503     {
504       const MEDCouplingUMesh *elt((*it).getMesh());
505       if(elt)
506         {
507           int myMeshDim(elt->getMeshDimension());
508           meshDim=std::max(meshDim,myMeshDim);
509           ms2[myMeshDim].push_back(*it);
510         }
511     }
512   if(ms2.empty())
513     return ;
514   for(std::map<int, std::vector< MicroField > >::const_iterator it=ms2.begin();it!=ms2.end();it++)
515     {
516       const std::vector< MicroField >& vs((*it).second);
517       if(vs.size()==1)
518         {
519           PutAtLevelDealOrder(mfd,(*it).first-meshDim,vs[0],timeStep,tsId);
520         }
521       else
522         {
523           MicroField merge(vs);
524           PutAtLevelDealOrder(mfd,(*it).first-meshDim,merge,timeStep,tsId);
525         }
526     }
527 }
528
529 DataArrayDouble *BuildCoordsFrom(vtkPointSet *ds)
530 {
531   if(!ds)
532     throw MZCException("BuildCoordsFrom : internal error !");
533   vtkPoints *pts(ds->GetPoints());
534   if(!pts)
535     throw MZCException("BuildCoordsFrom : internal error 2 !");
536   vtkDataArray *data(pts->GetData());
537   if(!data)
538     throw MZCException("BuildCoordsFrom : internal error 3 !");
539   return ConvertVTKArrayToMCArrayDoubleForced(data);
540 }
541
542 void AddNodeFields(MEDFileData *mfd, vtkDataSetAttributes *dsa, double timeStep, int tsId)
543 {
544   if(!mfd || !dsa)
545     throw MZCException("AddNodeFields : internal error !");
546   MEDFileMesh *mm(mfd->getMeshes()->getMeshAtPos(0));
547   MEDFileUMesh *mmu(dynamic_cast<MEDFileUMesh *>(mm));
548   if(!mmu)
549     throw MZCException("AddNodeFields : internal error 2 !");
550   MCAuto<MEDCouplingUMesh> mesh;
551   if(!mmu->getNonEmptyLevels().empty())
552     mesh=mmu->getMeshAtLevel(0);
553   else
554     {
555       mesh=MEDCouplingUMesh::Build0DMeshFromCoords(mmu->getCoords());
556       mesh->setName(mmu->getName());
557     }
558   int nba(dsa->GetNumberOfArrays());
559   for(int i=0;i<nba;i++)
560     {
561       vtkDataArray *arr(dsa->GetArray(i));
562       const char *name(arr->GetName());
563       if(!arr)
564         continue;
565       MCAuto<DataArray> da(ConvertVTKArrayToMCArray(arr));
566       da->setName(name);
567       AppendMCFieldFrom(MEDCoupling::ON_NODES,mesh,mfd,da,NULL,timeStep,tsId);
568     }
569 }
570
571 std::vector<MCAuto<DataArray> > AddPartFields(const DataArrayIdType *part, vtkDataSetAttributes *dsa)
572 {
573   std::vector< MCAuto<DataArray> > ret;
574   if(!dsa)
575     return ret;
576   int nba(dsa->GetNumberOfArrays());
577   for(int i=0;i<nba;i++)
578     {
579       vtkDataArray *arr(dsa->GetArray(i));
580       if(!arr)
581         continue;
582       const char *name(arr->GetName());
583       //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
584       //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
585       MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
586       if(part)
587         mcarr=mcarr->selectByTupleId(part->begin(),part->end());
588       mcarr->setName(name);
589       ret.push_back(mcarr);
590     }
591   return ret;
592 }
593
594 std::vector<MCAuto<DataArray> > AddPartFields2(int bg, int end, vtkDataSetAttributes *dsa)
595 {
596   std::vector< MCAuto<DataArray> > ret;
597   if(!dsa)
598     return ret;
599   int nba(dsa->GetNumberOfArrays());
600   for(int i=0;i<nba;i++)
601     {
602       vtkDataArray *arr(dsa->GetArray(i));
603       if(!arr)
604         continue;
605       const char *name(arr->GetName());
606       //int nbCompo(arr->GetNumberOfComponents()); // todo: unused
607       //vtkIdType nbTuples(arr->GetNumberOfTuples()); // todo: unused
608       MCAuto<DataArray> mcarr(ConvertVTKArrayToMCArray(arr));
609       mcarr=mcarr->selectByTupleIdSafeSlice(bg,end,1);
610       mcarr->setName(name);
611       ret.push_back(mcarr);
612     }
613   return ret;
614 }
615
616 void ConvertFromRectilinearGrid(MEDFileData *ret, vtkRectilinearGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
617 {
618   if(!ds || !ret)
619     throw MZCException("ConvertFromRectilinearGrid : internal error !");
620   //
621   MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
622   ret->setMeshes(meshes);
623   MCAuto<MEDFileFields> fields(MEDFileFields::New());
624   ret->setFields(fields);
625   //
626   MCAuto<MEDFileCMesh> cmesh(MEDFileCMesh::New());
627   meshes->pushMesh(cmesh);
628   MCAuto<MEDCouplingCMesh> cmeshmc(MEDCouplingCMesh::New());
629   vtkDataArray *cx(ds->GetXCoordinates()),*cy(ds->GetYCoordinates()),*cz(ds->GetZCoordinates());
630   if(cx)
631     {
632       MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cx));
633       cmeshmc->setCoordsAt(0,arr);
634     }
635   if(cy)
636     {
637       MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cy));
638       cmeshmc->setCoordsAt(1,arr);
639     }
640   if(cz)
641     {
642       MCAuto<DataArrayDouble> arr(ConvertVTKArrayToMCArrayDoubleForced(cz));
643       cmeshmc->setCoordsAt(2,arr);
644     }
645   std::string meshName(GetMeshNameWithContext(context));
646   cmeshmc->setName(meshName);
647   cmesh->setMesh(cmeshmc);
648   std::vector<MCAuto<DataArray> > cellFs(AddPartFields(0,ds->GetCellData()));
649   for(std::vector<MCAuto<DataArray> >::const_iterator it=cellFs.begin();it!=cellFs.end();it++)
650     {
651       MCAuto<DataArray> da(*it);
652       AppendMCFieldFrom(MEDCoupling::ON_CELLS,cmeshmc,ret,da,NULL,timeStep,tsId);
653     }
654   std::vector<MCAuto<DataArray> > nodeFs(AddPartFields(0,ds->GetPointData()));
655   for(std::vector<MCAuto<DataArray> >::const_iterator it=nodeFs.begin();it!=nodeFs.end();it++)
656     {
657       MCAuto<DataArray> da(*it);
658       AppendMCFieldFrom(MEDCoupling::ON_NODES,cmeshmc,ret,da,NULL,timeStep,tsId);
659     }
660 }
661
662 void ConvertFromPolyData(MEDFileData *ret, vtkPolyData *ds, const std::vector<int>& context, double timeStep, int tsId)
663 {
664   if(!ds || !ret)
665     throw MZCException("ConvertFromPolyData : internal error !");
666   //
667   MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
668   ret->setMeshes(meshes);
669   MCAuto<MEDFileFields> fields(MEDFileFields::New());
670   ret->setFields(fields);
671   //
672   MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
673   meshes->pushMesh(umesh);
674   MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
675   umesh->setCoords(coords);
676   umesh->setName(GetMeshNameWithContext(context));
677   //
678   int offset(0);
679   std::vector< MicroField > ms;
680   vtkCellArray *cd(ds->GetVerts());
681   if(cd)
682     {
683       MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cd,coords,0,INTERP_KERNEL::NORM_POINT1));
684       if((const MEDCouplingUMesh *)subMesh)
685         {
686           std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
687           offset+=subMesh->getNumberOfCells();
688           ms.push_back(MicroField(subMesh,cellFs));
689         }
690     }
691   vtkCellArray *cc(ds->GetLines());
692   if(cc)
693     {
694       MCAuto<MEDCouplingUMesh> subMesh;
695       try
696         {
697           subMesh=BuildMeshFromCellArray(cc,coords,1,INTERP_KERNEL::NORM_SEG2);
698         }
699       catch(INTERP_KERNEL::Exception& e)
700         {
701           std::ostringstream oss; oss << "MEDWriter does not manage polyline cell type because MED file format does not support it ! Maybe it is the source of the problem ? The cause of this exception was " << e.what() << std::endl;
702           throw INTERP_KERNEL::Exception(oss.str());
703         }
704       if((const MEDCouplingUMesh *)subMesh)
705         {
706           std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
707           offset+=subMesh->getNumberOfCells();
708           ms.push_back(MicroField(subMesh,cellFs));
709         }
710     }
711   vtkCellArray *cb(ds->GetPolys());
712   if(cb)
713     {
714       MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArray(cb,coords,2,INTERP_KERNEL::NORM_POLYGON));
715       if((const MEDCouplingUMesh *)subMesh)
716         {
717           std::vector<MCAuto<DataArray> > cellFs(AddPartFields2(offset,offset+subMesh->getNumberOfCells(),ds->GetCellData()));
718           offset+=subMesh->getNumberOfCells();
719           ms.push_back(MicroField(subMesh,cellFs));
720         }
721     }
722   vtkCellArray *ca(ds->GetStrips());
723   if(ca)
724     {
725       MCAuto<DataArrayIdType> ids;
726       MCAuto<MEDCouplingUMesh> subMesh(BuildMeshFromCellArrayTriangleStrip(ca,coords,ids));
727       if((const MEDCouplingUMesh *)subMesh)
728         {
729           std::vector<MCAuto<DataArray> > cellFs(AddPartFields(ids,ds->GetCellData()));
730           offset+=subMesh->getNumberOfCells();
731           ms.push_back(MicroField(subMesh,cellFs));
732         }
733     }
734   AssignSingleGTMeshes(ret,ms,timeStep,tsId);
735   AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
736 }
737
738 void ConvertFromUnstructuredGrid(MEDFileData *ret, vtkUnstructuredGrid *ds, const std::vector<int>& context, double timeStep, int tsId)
739 {
740   if(!ds || !ret)
741     throw MZCException("ConvertFromUnstructuredGrid : internal error !");
742   //
743   MCAuto<MEDFileMeshes> meshes(MEDFileMeshes::New());
744   ret->setMeshes(meshes);
745   MCAuto<MEDFileFields> fields(MEDFileFields::New());
746   ret->setFields(fields);
747   //
748   MCAuto<MEDFileUMesh> umesh(MEDFileUMesh::New());
749   meshes->pushMesh(umesh);
750   MCAuto<DataArrayDouble> coords(BuildCoordsFrom(ds));
751   umesh->setCoords(coords);
752   umesh->setName(GetMeshNameWithContext(context));
753   vtkIdType nbCells(ds->GetNumberOfCells());
754   vtkCellArray *ca(ds->GetCells());
755   if(!ca)
756     return ;
757   //vtkIdType nbEnt(ca->GetNumberOfConnectivityEntries()); // todo: unused
758   //vtkIdType *caPtr(ca->GetData()->GetPointer(0)); // todo: unused
759   vtkUnsignedCharArray *ct(ds->GetCellTypesArray());
760   if(!ct)
761     throw MZCException("ConvertFromUnstructuredGrid : internal error");
762   vtkIdTypeArray *cla(ds->GetCellLocationsArray());
763   //const vtkIdType *claPtr(cla->GetPointer(0)); // todo: unused
764   if(!cla)
765     throw MZCException("ConvertFromUnstructuredGrid : internal error 2");
766   const unsigned char *ctPtr(ct->GetPointer(0));
767   std::map<int,int> m(ComputeMapOfType());
768   MCAuto<DataArrayInt> lev(DataArrayInt::New()) ;  lev->alloc(nbCells,1);
769   int *levPtr(lev->getPointer());
770   for(vtkIdType i=0;i<nbCells;i++)
771     {
772       std::map<int,int>::iterator it(m.find(ctPtr[i]));
773       if(it!=m.end())
774         {
775           const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)(*it).second));
776           levPtr[i]=cm.getDimension();
777         }
778       else
779         {
780           if(ctPtr[i]==VTK_POLY_VERTEX)
781             {
782               const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(INTERP_KERNEL::NORM_POINT1));
783               levPtr[i]=cm.getDimension();
784             }
785           else
786             {
787               std::ostringstream oss; oss << "ConvertFromUnstructuredGrid : at pos #" << i << " unrecognized VTK cell with type =" << ctPtr[i];
788               throw MZCException(oss.str());
789             }
790         }
791     }
792   MCAuto<DataArrayInt> levs(lev->getDifferentValues());
793   std::vector< MicroField > ms;
794   //vtkIdTypeArray *faces(ds->GetFaces()),*faceLoc(ds->GetFaceLocations()); // todo: unused
795   for(const int *curLev=levs->begin();curLev!=levs->end();curLev++)
796     {
797       MCAuto<MEDCouplingUMesh> m0(MEDCouplingUMesh::New("",*curLev));
798       m0->setCoords(coords); m0->allocateCells();
799       MCAuto<DataArrayIdType> cellIdsCurLev(lev->findIdsEqual(*curLev));
800       for(const mcIdType *cellId=cellIdsCurLev->begin();cellId!=cellIdsCurLev->end();cellId++)
801         {
802           int vtkType(ds->GetCellType(*cellId));
803           std::map<int,int>::iterator it(m.find(vtkType));
804           INTERP_KERNEL::NormalizedCellType ct=it!=m.end()?(INTERP_KERNEL::NormalizedCellType)((*it).second):INTERP_KERNEL::NORM_POINT1;
805           if(ct!=INTERP_KERNEL::NORM_POLYHED && vtkType!=VTK_POLY_VERTEX)
806             {
807               vtkIdType sz(0);
808               const vtkIdType *pts(nullptr);
809               ds->GetCellPoints(*cellId, sz, pts);
810               std::vector<mcIdType> conn2(pts,pts+sz);
811               m0->insertNextCell(ct,sz,conn2.data());
812             }
813           else if(ct==INTERP_KERNEL::NORM_POLYHED)
814             {
815               // # de faces du polyèdre
816               vtkIdType nbOfFaces(0);
817               // connectivé des faces (numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...)
818               const vtkIdType *facPtr(nullptr);
819               ds->GetFaceStream(*cellId, nbOfFaces, facPtr);
820               std::vector<mcIdType> conn;
821               for(vtkIdType k=0;k<nbOfFaces;k++)
822                 {
823                   vtkIdType nbOfNodesInFace(*facPtr++);
824                   std::copy(facPtr,facPtr+nbOfNodesInFace,std::back_inserter(conn));
825                   if(k<nbOfFaces-1)
826                     conn.push_back(-1);
827                   facPtr+=nbOfNodesInFace;
828                 }
829               m0->insertNextCell(ct,ToIdType(conn.size()),&conn[0]);
830             }
831           else
832             {
833               vtkIdType sz(0);
834               const vtkIdType *pts(nullptr);
835               ds->GetCellPoints(*cellId, sz, pts);
836               if(sz!=1)
837                 throw MZCException("ConvertFromUnstructuredGrid : non single poly vertex not managed by MED !");
838               m0->insertNextCell(ct,1,(const mcIdType*)pts);
839             }
840         }
841       std::vector<MCAuto<DataArray> > cellFs(AddPartFields(cellIdsCurLev,ds->GetCellData()));
842       ms.push_back(MicroField(m0,cellFs));
843     }
844   AssignSingleGTMeshes(ret,ms,timeStep,tsId);
845   AddNodeFields(ret,ds->GetPointData(),timeStep,tsId);
846 }
847
848 ///////////////////
849
850 void WriteMEDFileFromVTKDataSet(MEDFileData *mfd, vtkDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
851 {
852   if(!ds || !mfd)
853     throw MZCException("Internal error in WriteMEDFileFromVTKDataSet.");
854   vtkPolyData *ds2(vtkPolyData::SafeDownCast(ds));
855   vtkUnstructuredGrid *ds3(vtkUnstructuredGrid::SafeDownCast(ds));
856   vtkRectilinearGrid *ds4(vtkRectilinearGrid::SafeDownCast(ds));
857   if(ds2)
858     {
859       ConvertFromPolyData(mfd,ds2,context,timeStep,tsId);
860     }
861   else if(ds3)
862     {
863       ConvertFromUnstructuredGrid(mfd,ds3,context,timeStep,tsId);
864     }
865   else if(ds4)
866     {
867       ConvertFromRectilinearGrid(mfd,ds4,context,timeStep,tsId);
868     }
869   else
870     throw MZCException("Unrecognized vtkDataSet ! Sorry ! Try to convert it to UnstructuredGrid to be able to write it !");
871 }
872
873 void WriteMEDFileFromVTKMultiBlock(MEDFileData *mfd, vtkMultiBlockDataSet *ds, const std::vector<int>& context, double timeStep, int tsId)
874 {
875   if(!ds || !mfd)
876     throw MZCException("Internal error in WriteMEDFileFromVTKMultiBlock.");
877   int nbBlocks(ds->GetNumberOfBlocks());
878   if(nbBlocks==1 && context.empty())
879     {
880       vtkDataObject *uniqueElt(ds->GetBlock(0));
881       if(!uniqueElt)
882         throw MZCException("Unique elt in multiblock is NULL !");
883       vtkDataSet *uniqueEltc(vtkDataSet::SafeDownCast(uniqueElt));
884       if(uniqueEltc)
885         {
886           WriteMEDFileFromVTKDataSet(mfd,uniqueEltc,context,timeStep,tsId);
887           return ;
888         }
889     }
890   for(int i=0;i<nbBlocks;i++)
891     {
892       vtkDataObject *elt(ds->GetBlock(i));
893       std::vector<int> context2;
894       context2.push_back(i);
895       if(!elt)
896         {
897           std::ostringstream oss; oss << "In context ";
898           std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
899           oss << " at pos #" << i << " elt is NULL !";
900           throw MZCException(oss.str());
901         }
902       vtkDataSet *elt1(vtkDataSet::SafeDownCast(elt));
903       if(elt1)
904         {
905           WriteMEDFileFromVTKDataSet(mfd,elt1,context,timeStep,tsId);
906           continue;
907         }
908       vtkMultiBlockDataSet *elt2(vtkMultiBlockDataSet::SafeDownCast(elt));
909       if(elt2)
910         {
911           WriteMEDFileFromVTKMultiBlock(mfd,elt2,context,timeStep,tsId);
912           continue;
913         }
914       std::ostringstream oss; oss << "In context ";
915       std::copy(context.begin(),context.end(),std::ostream_iterator<int>(oss," "));
916       oss << " at pos #" << i << " elt not recognized data type !";
917       throw MZCException(oss.str());
918     }
919 }
920
921 void WriteMEDFileFromVTKGDS(MEDFileData *mfd, vtkDataObject *input, double timeStep, int tsId)
922 {
923   if(!input || !mfd)
924     throw MZCException("WriteMEDFileFromVTKGDS : internal error !");
925   std::vector<int> context;
926   vtkDataSet *input1(vtkDataSet::SafeDownCast(input));
927   if(input1)
928     {
929       WriteMEDFileFromVTKDataSet(mfd,input1,context,timeStep,tsId);
930       return ;
931     }
932   vtkMultiBlockDataSet *input2(vtkMultiBlockDataSet::SafeDownCast(input));
933   if(input2)
934     {
935       WriteMEDFileFromVTKMultiBlock(mfd,input2,context,timeStep,tsId);
936       return ;
937     }
938   throw MZCException("WriteMEDFileFromVTKGDS : not recognized data type !");
939 }
940
941 void PutFamGrpInfoIfAny(MEDFileData *mfd, const std::string& meshName, const std::vector<Grp>& groups, const std::vector<Fam>& fams)
942 {
943   if(!mfd)
944     return ;
945   if(meshName.empty())
946     return ;
947   MEDFileMeshes *meshes(mfd->getMeshes());
948   if(!meshes)
949     return ;
950   if(meshes->getNumberOfMeshes()!=1)
951     return ;
952   MEDFileMesh *mm(meshes->getMeshAtPos(0));
953   if(!mm)
954     return ;
955   mm->setName(meshName);
956   for(std::vector<Fam>::const_iterator it=fams.begin();it!=fams.end();it++)
957     mm->setFamilyId((*it).getName(),(*it).getID());
958   for(std::vector<Grp>::const_iterator it=groups.begin();it!=groups.end();it++)
959     mm->setFamiliesOnGroup((*it).getName(),(*it).getFamilies());
960   MEDFileFields *fields(mfd->getFields());
961   if(!fields)
962     return ;
963   for(int i=0;i<fields->getNumberOfFields();i++)
964     {
965       MEDFileAnyTypeFieldMultiTS *fmts(fields->getFieldAtPos(i));
966       if(!fmts)
967         continue;
968       fmts->setMeshName(meshName);
969     }
970 }