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