Salome HOME
PR: debug polygons and merge nodes
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
1 #define CHRONODEF
2 #include "SMDS_UnstructuredGrid.hxx"
3 #include "SMDS_Mesh.hxx"
4 #include "SMDS_MeshInfo.hxx"
5 #include "SMDS_Downward.hxx"
6
7 #include "utilities.h"
8
9 #include <vtkCellArray.h>
10 #include <vtkCellLinks.h>
11 #include <vtkIdTypeArray.h>
12 #include <vtkUnsignedCharArray.h>
13
14 #include <list>
15
16 using namespace std;
17
18 SMDS_CellLinks* SMDS_CellLinks::New()
19 {
20   MESSAGE("SMDS_CellLinks::New");
21   return new SMDS_CellLinks();
22 }
23
24 vtkCellLinks::Link* SMDS_CellLinks::ResizeL(vtkIdType sz)
25 {
26   return vtkCellLinks::Resize(sz);
27 }
28
29 vtkIdType SMDS_CellLinks::GetLinksSize()
30 {
31   return this->Size;
32 }
33
34 SMDS_CellLinks::SMDS_CellLinks() :
35   vtkCellLinks()
36 {
37 }
38
39 SMDS_CellLinks::~SMDS_CellLinks()
40 {
41 }
42
43 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
44 {
45   MESSAGE("SMDS_UnstructuredGrid::New");
46   return new SMDS_UnstructuredGrid();
47 }
48
49 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
50   vtkUnstructuredGrid()
51 {
52   _cellIdToDownId.clear();
53   _downTypes.clear();
54   _downArray.clear();
55   _mesh = 0;
56 }
57
58 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
59 {
60 }
61
62 unsigned long SMDS_UnstructuredGrid::GetMTime()
63 {
64   unsigned long mtime = vtkUnstructuredGrid::GetMTime();
65   MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
66   return mtime;
67 }
68
69 void SMDS_UnstructuredGrid::Update()
70 {
71   MESSAGE("SMDS_UnstructuredGrid::Update");
72   return vtkUnstructuredGrid::Update();
73 }
74
75 void SMDS_UnstructuredGrid::UpdateInformation()
76 {
77   MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
78   return vtkUnstructuredGrid::UpdateInformation();
79 }
80
81 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
82 {
83   // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
84   //MESSAGE("*********************** SMDS_UnstructuredGrid::GetPoints " << this->Points << " " << vtkUnstructuredGrid::GetPoints());
85   return this->Points;
86 }
87
88 //#ifdef VTK_HAVE_POLYHEDRON
89 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
90 {
91   if (type != VTK_POLYHEDRON)
92     return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
93
94   // --- type = VTK_POLYHEDRON
95   //MESSAGE("InsertNextLinkedCell VTK_POLYHEDRON");
96   int cellid = this->InsertNextCell(type, npts, pts);
97
98   set<vtkIdType> setOfNodes;
99   setOfNodes.clear();
100   int nbfaces = npts;
101   int i = 0;
102   for (int nf = 0; nf < nbfaces; nf++)
103     {
104       int nbnodes = pts[i];
105       i++;
106       for (int k = 0; k < nbnodes; k++)
107         {
108           //MESSAGE(" cell " << cellid << " face " << nf << " node " << pts[i]);
109           setOfNodes.insert(pts[i]);
110           i++;
111         }
112     }
113
114   set<vtkIdType>::iterator it = setOfNodes.begin();
115   for (; it != setOfNodes.end(); ++it)
116     {
117       //MESSAGE("reverse link for node " << *it << " cell " << cellid);
118       this->Links->ResizeCellList(*it, 1);
119       this->Links->AddCellReference(cellid, *it);
120     }
121
122   return cellid;
123 }
124 //#endif
125
126 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
127 {
128   _mesh = mesh;
129 }
130
131 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
132                                         std::vector<int>& idCellsOldToNew, int newCellSize)
133 {
134   // TODO utiliser mieux vtk pour faire plus simple (plus couteux ?)
135
136   MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
137   int startHole = 0;
138   int endHole = 0;
139   int startBloc = 0;
140   int endBloc = 0;
141   int alreadyCopied = 0;
142   int holes = 0;
143
144   typedef enum
145   {
146     lookHoleStart, lookHoleEnd, lookBlocEnd
147   } enumState;
148   enumState compactState = lookHoleStart;
149
150   //    if (this->Links)
151   //    {
152   //            this->Links->UnRegister(this);
153   //            this->Links = 0;
154   //    }
155
156   // --- if newNodeSize, create a new compacted vtkPoints
157
158   vtkPoints *newPoints = 0;
159   if (newNodeSize)
160     {
161       MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
162       newPoints = vtkPoints::New();
163       newPoints->Initialize();
164       newPoints->Allocate(newNodeSize);
165       newPoints->SetNumberOfPoints(newNodeSize);
166       int oldNodeSize = idNodesOldToNew.size();
167
168       for (int i = 0; i < oldNodeSize; i++)
169         {
170           //MESSAGE("                                    " << i << " " << idNodesOldToNew[i]);
171           switch (compactState)
172           {
173             case lookHoleStart:
174               if (idNodesOldToNew[i] < 0)
175                 {
176                   MESSAGE("-------------- newNodeSize, startHole " << i << " " << oldNodeSize);
177                   startHole = i;
178                   if (!alreadyCopied) // copy the first bloc
179                     {
180                       MESSAGE("--------- copy first nodes before hole " << i << " " << oldNodeSize);
181                       copyNodes(newPoints, idNodesOldToNew, alreadyCopied, 0, startHole);
182                     }
183                   compactState = lookHoleEnd;
184                 }
185               break;
186             case lookHoleEnd:
187               if (idNodesOldToNew[i] >= 0)
188                 {
189                   MESSAGE("-------------- newNodeSize, endHole " << i << " " << oldNodeSize);
190                   endHole = i;
191                   startBloc = i;
192                   compactState = lookBlocEnd;
193                 }
194               break;
195             case lookBlocEnd:
196               if (idNodesOldToNew[i] < 0)
197                 endBloc = i; // see nbPoints below
198               else if (i == (oldNodeSize - 1))
199                 endBloc = i + 1;
200               if (endBloc)
201                 {
202                   MESSAGE("-------------- newNodeSize, endbloc " << endBloc << " " << oldNodeSize);
203                   copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
204                   compactState = lookHoleEnd;
205                   startHole = i;
206                   endHole = 0;
207                   startBloc = 0;
208                   endBloc = 0;
209                 }
210               break;
211           }
212         }
213       if (!alreadyCopied) // no hole, but shorter, no need to modify idNodesOldToNew
214         {
215           MESSAGE("------------- newNodeSize, shorter " << oldNodeSize);
216           copyNodes(newPoints, idNodesOldToNew, alreadyCopied, 0, newNodeSize);
217         }
218       newPoints->Squeeze();
219     }
220
221   // --- create new compacted Connectivity, Locations and Types
222
223   int oldCellSize = this->Types->GetNumberOfTuples();
224
225   vtkCellArray *newConnectivity = vtkCellArray::New();
226   newConnectivity->Initialize();
227   int oldCellDataSize = this->Connectivity->GetData()->GetSize();
228   newConnectivity->Allocate(oldCellDataSize);
229   MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
230
231   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
232   newTypes->Initialize();
233   newTypes->SetNumberOfValues(newCellSize);
234
235   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
236   newLocations->Initialize();
237   newLocations->SetNumberOfValues(newCellSize);
238
239   startHole = 0;
240   endHole = 0;
241   startBloc = 0;
242   endBloc = 0;
243   alreadyCopied = 0;
244   holes = 0;
245   compactState = lookHoleStart;
246
247   // TODO some polyhedron may be huge (only in some tests)
248   vtkIdType tmpid[NBMAXNODESINCELL];
249   vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
250
251   for (int i = 0; i < oldCellSize; i++)
252     {
253       switch (compactState)
254       {
255         case lookHoleStart:
256           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
257             {
258               MESSAGE(" -------- newCellSize, startHole " << i << " " << oldCellSize);
259               startHole = i;
260               compactState = lookHoleEnd;
261               if (!alreadyCopied) // copy the first bloc
262                 {
263                   MESSAGE("--------- copy first bloc before hole " << i << " " << oldCellSize);
264                   copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
265                            alreadyCopied, 0, startHole);
266                 }
267             }
268           break;
269         case lookHoleEnd:
270           if (this->Types->GetValue(i) != VTK_EMPTY_CELL)
271             {
272               MESSAGE(" -------- newCellSize, EndHole " << i << " " << oldCellSize);
273               endHole = i;
274               startBloc = i;
275               compactState = lookBlocEnd;
276               holes += endHole - startHole;
277             }
278           break;
279         case lookBlocEnd:
280           endBloc = 0;
281           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
282             endBloc = i;
283           else if (i == (oldCellSize - 1))
284             endBloc = i + 1;
285           if (endBloc)
286             {
287               MESSAGE(" -------- newCellSize, endBloc " << endBloc << " " << oldCellSize);
288               copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
289                        alreadyCopied, startBloc, endBloc);
290               compactState = lookHoleEnd;
291             }
292           break;
293       }
294     }
295   if (!alreadyCopied) // no hole, but shorter
296     {
297       MESSAGE(" -------- newCellSize, shorter " << oldCellSize);
298       copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell, alreadyCopied, 0,
299                oldCellSize);
300     }
301
302   newConnectivity->Squeeze();
303   //newTypes->Squeeze();
304   //newLocations->Squeeze();
305
306   if (newNodeSize)
307     {
308       MESSAGE("------- newNodeSize, setPoints");
309       this->SetPoints(newPoints);
310       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
311     }
312
313   if (this->FaceLocations)
314     {
315       vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
316       newFaceLocations->Initialize();
317       newFaceLocations->Allocate(newTypes->GetSize());
318       vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
319       newFaces->Initialize();
320       newFaces->Allocate(this->Faces->GetSize());
321       for (int i = 0; i < oldCellSize; i++)
322         {
323           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
324             continue;
325           int newCellId = idCellsOldToNew[i];
326           if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
327             {
328                newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
329                int oldFaceLoc = this->FaceLocations->GetValue(i);
330                int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
331                newFaces->InsertNextValue(nCellFaces);
332                for (int n=0; n<nCellFaces; n++)
333                  {
334                    int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
335                    newFaces->InsertNextValue(nptsInFace);
336                    for (int k=0; k<nptsInFace; k++)
337                      {
338                        int oldpt = this->Faces->GetValue(oldFaceLoc++);
339                        newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
340                      }
341                  }
342             }
343           else
344             {
345                newFaceLocations->InsertNextValue(-1);
346             }
347         }
348       newFaceLocations->Squeeze();
349       newFaces->Squeeze();
350       newFaceLocations->Register(this);
351       newFaces->Register(this);
352       this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
353       newFaceLocations->Delete();
354       newFaces->Delete();
355     }
356   else
357     this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
358
359   newTypes->Delete();
360   newLocations->Delete();
361   newConnectivity->Delete();
362   this->BuildLinks();
363 }
364
365 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
366                                       int start, int end)
367 {
368   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
369   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
370   void *source = this->Points->GetVoidPointer(3 * start);
371   int nbPoints = end - start;
372   if (nbPoints > 0)
373     {
374       memcpy(target, source, 3 * sizeof(float) * nbPoints);
375       for (int j = start; j < end; j++)
376         idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
377         //idNodesOldToNew[alreadyCopied++] = idNodesOldToNew[j]; // new vtkId --> old SMDS id
378     }
379 }
380
381 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
382                                      std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
383                                      vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
384                                      int start, int end)
385 {
386   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
387   for (int j = start; j < end; j++)
388     {
389       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
390       idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
391       vtkIdType oldLoc = this->Locations->GetValue(j);
392       vtkIdType nbpts;
393       vtkIdType *oldPtsCell = 0;
394       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
395       assert(nbpts < NBMAXNODESINCELL);
396       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
397       for (int l = 0; l < nbpts; l++)
398         {
399           int oldval = oldPtsCell[l];
400           pointsCell[l] = idNodesOldToNew[oldval];
401           //MESSAGE("   " << oldval << " " << pointsCell[l]);
402         }
403       int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
404       int newLoc = newConnectivity->GetInsertLocation(nbpts);
405       //MESSAGE(newcnt << " " << newLoc);
406       newLocations->SetValue(alreadyCopied, newLoc);
407       alreadyCopied++;
408     }
409 }
410
411 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
412 {
413   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
414   return _cellIdToDownId[vtkCellId];
415 }
416
417 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
418 {
419   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
420   _cellIdToDownId[vtkCellId] = downId;
421 }
422
423 /*! Build downward connectivity: to do only when needed because heavy memory load.
424  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
425  *
426  */
427 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
428 {
429   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
430   // TODO calcul partiel sans edges
431
432   // --- erase previous data if any
433
434   for (int i = 0; i < _downArray.size(); i++)
435     {
436       if (_downArray[i])
437         delete _downArray[i];
438       _downArray[i] = 0;
439     }
440   _cellIdToDownId.clear();
441
442   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
443
444   _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
445
446   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
447   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
448   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
449   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
450   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
451   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
452   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
453   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
454   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
455   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
456   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
457   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
458   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
459   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
460
461   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
462
463   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
464
465   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
466   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
467   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
468   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
469   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
470   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
471   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
472   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
473
474   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
475   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
476   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
477   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
478   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
479   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
480
481   int GuessSize[VTK_QUADRATIC_TETRA];
482   GuessSize[VTK_LINE] = nbLineGuess;
483   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
484   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
485   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
486   GuessSize[VTK_QUAD] = nbLinQuadGuess;
487   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
488   GuessSize[VTK_TETRA] = nbLinTetra;
489   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
490   GuessSize[VTK_PYRAMID] = nbLinPyra;
491   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
492   GuessSize[VTK_WEDGE] = nbLinPrism;
493   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
494   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
495   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
496
497   _downArray[VTK_LINE]->allocate(nbLineGuess);
498   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
499   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
500   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
501   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
502   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
503   _downArray[VTK_TETRA]->allocate(nbLinTetra);
504   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
505   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
506   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
507   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
508   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
509   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
510   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
511
512   // --- iteration on vtkUnstructuredGrid cells, only faces
513   //     for each vtk face:
514   //       create a downward face entry with its downward id.
515   //       compute vtk volumes, create downward volumes entry.
516   //       mark face in downward volumes
517   //       mark volumes in downward face
518
519   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
520   int cellSize = this->Types->GetNumberOfTuples();
521   _cellIdToDownId.resize(cellSize, -1);
522
523   for (int i = 0; i < cellSize; i++)
524     {
525       int vtkFaceType = this->GetCellType(i);
526       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
527         {
528           int vtkFaceId = i;
529           //ASSERT(_downArray[vtkFaceType]);
530           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
531           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
532           downFace->setTempNodes(connFaceId, vtkFaceId);
533           int vols[2] = { -1, -1 };
534           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
535           //MESSAGE("nbVolumes="<< nbVolumes);
536           for (int ivol = 0; ivol < nbVolumes; ivol++)
537             {
538               int vtkVolId = vols[ivol];
539               int vtkVolType = this->GetCellType(vtkVolId);
540               //ASSERT(_downArray[vtkVolType]);
541               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
542               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
543               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
544               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
545             }
546         }
547     }
548
549   // --- iteration on vtkUnstructuredGrid cells, only volumes
550   //     for each vtk volume:
551   //       create downward volumes entry if not already done
552   //       build a temporary list of faces described with their nodes
553   //       for each face
554   //         compute the vtk volumes containing this face
555   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
556   //         if not, create a downward face entry (resizing of structure required)
557   //         (the downward faces store a temporary list of nodes to ease the comparison)
558   //         create downward volumes entry if not already done
559   //         mark volumes in downward face
560   //         mark face in downward volumes
561
562   CHRONOSTOP(20);
563   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
564
565   for (int i = 0; i < cellSize; i++)
566     {
567       int vtkType = this->GetCellType(i);
568       if (SMDS_Downward::getCellDimension(vtkType) == 3)
569         {
570           //CHRONO(31);
571           int vtkVolId = i;
572           // MESSAGE("vtk volume " << vtkVolId);
573           //ASSERT(_downArray[vtkType]);
574           int connVolId = _downArray[vtkType]->addCell(vtkVolId);
575
576           // --- find all the faces of the volume, describe the faces by their nodes
577
578           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
579           ListElemByNodesType facesWithNodes;
580           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
581           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
582           //CHRONOSTOP(31);
583           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
584             {
585               // --- find the volumes containing the face
586
587               //CHRONO(32);
588               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
589               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
590               int vols[2] = { -1, -1 };
591               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
592               int lg = facesWithNodes.elems[iface].nbNodes;
593               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
594               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
595
596               // --- check if face is registered in the volumes
597               //CHRONOSTOP(32);
598
599               //CHRONO(33);
600               int connFaceId = -1;
601               for (int ivol = 0; ivol < nbVolumes; ivol++)
602                 {
603                   int vtkVolId2 = vols[ivol];
604                   int vtkVolType = this->GetCellType(vtkVolId2);
605                   //ASSERT(_downArray[vtkVolType]);
606                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
607                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
608                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
609                   if (connFaceId >= 0)
610                     break; // --- face already created
611                 }//CHRONOSTOP(33);
612
613               // --- if face is not registered in the volumes, create face
614
615               //CHRONO(34);
616               if (connFaceId < 0)
617                 {
618                   connFaceId = _downArray[vtkFaceType]->addCell();
619                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
620                 }//CHRONOSTOP(34);
621
622               // --- mark volumes in downward face and mark face in downward volumes
623
624               //CHRONO(35);
625               for (int ivol = 0; ivol < nbVolumes; ivol++)
626                 {
627                   int vtkVolId2 = vols[ivol];
628                   int vtkVolType = this->GetCellType(vtkVolId2);
629                   //ASSERT(_downArray[vtkVolType]);
630                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
631                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
632                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
633                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
634                 }//CHRONOSTOP(35);
635             }
636         }
637     }
638
639   // --- iteration on vtkUnstructuredGrid cells, only edges
640   //     for each vtk edge:
641   //       create downward edge entry
642   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
643   //       find downward faces
644   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
645   //       mark edge in downward faces
646   //       mark faces in downward edge
647
648   CHRONOSTOP(21);
649   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
650
651   for (int i = 0; i < cellSize; i++)
652     {
653       int vtkEdgeType = this->GetCellType(i);
654       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
655         {
656           int vtkEdgeId = i;
657           //ASSERT(_downArray[vtkEdgeType]);
658           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
659           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
660           downEdge->setNodes(connEdgeId, vtkEdgeId);
661           vector<int> vtkIds;
662           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
663           int downFaces[1000];
664           unsigned char downTypes[1000];
665           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
666           for (int n = 0; n < nbDownFaces; n++)
667             {
668               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
669               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
670             }
671         }
672     }
673
674   // --- iteration on downward faces (they are all listed now)
675   //     for each downward face:
676   //       build a temporary list of edges with their ordered list of nodes
677   //       for each edge:
678   //         find all the vtk cells containing this edge
679   //         then identify all the downward faces containing the edge, from the vtk cells
680   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
681   //         if not, create a downward edge entry with the node id's
682   //         mark edge in downward faces
683   //         mark downward faces in edge (size of list unknown, to be allocated)
684
685   CHRONOSTOP(22);CHRONO(23);
686
687   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
688     {
689       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
690         continue;
691
692       // --- find all the edges of the face, describe the edges by their nodes
693
694       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
695       int maxId = downFace->getMaxId();
696       for (int faceId = 0; faceId < maxId; faceId++)
697         {
698           //CHRONO(40);
699           ListElemByNodesType edgesWithNodes;
700           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
701           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
702
703           //CHRONOSTOP(40);
704           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
705             {
706
707               // --- check if the edge is already registered by exploration of the faces
708
709               //CHRONO(41);
710               vector<int> vtkIds;
711               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
712               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
713               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
714               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
715               //CHRONOSTOP(41);CHRONO(42);
716               int downFaces[1000];
717               unsigned char downTypes[1000];
718               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
719               //CHRONOSTOP(42);
720
721               //CHRONO(43);
722               int connEdgeId = -1;
723               for (int idf = 0; idf < nbDownFaces; idf++)
724                 {
725                   int faceId2 = downFaces[idf];
726                   int faceType = downTypes[idf];
727                   //ASSERT(_downArray[faceType]);
728                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
729                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
730                   if (connEdgeId >= 0)
731                     break; // --- edge already created
732                 }//CHRONOSTOP(43);
733
734               // --- if edge is not registered in the faces, create edge
735
736               if (connEdgeId < 0)
737                 {
738                   //CHRONO(44);
739                   connEdgeId = _downArray[vtkEdgeType]->addCell();
740                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
741                   //CHRONOSTOP(44);
742                 }
743
744               // --- mark faces in downward edge and mark edge in downward faces
745
746               //CHRONO(45);
747               for (int idf = 0; idf < nbDownFaces; idf++)
748                 {
749                   int faceId2 = downFaces[idf];
750                   int faceType = downTypes[idf];
751                   //ASSERT(_downArray[faceType]);
752                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
753                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
754                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
755                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
756                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
757                 }//CHRONOSTOP(45);
758             }
759         }
760     }
761
762   CHRONOSTOP(23);CHRONO(24);
763
764   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
765   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
766
767   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
768     {
769       if (SMDS_Downward *down = _downArray[vtkType])
770         {
771           down->compactStorage();
772         }
773     }
774
775   // --- Statistics
776
777   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
778     {
779       if (SMDS_Downward *down = _downArray[vtkType])
780         {
781           if (down->getMaxId())
782             {
783               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
784                   << GuessSize[vtkType] << " real: " << down->getMaxId());
785             }
786         }
787     }CHRONOSTOP(24);CHRONOSTOP(2);
788   counters::stats();
789 }
790
791 /*! Get the neighbors of a cell.
792  * Only the neighbors having the dimension of the cell are taken into account
793  * (neighbors of a volume are the volumes sharing a face with this volume,
794  *  neighbors of a face are the faces sharing an edge with this face...).
795  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
796  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
797  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
798  * @param vtkId the vtk id of the cell
799  * @return number of neighbors
800  */
801 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
802 {
803   int vtkType = this->GetCellType(vtkId);
804   int cellDim = SMDS_Downward::getCellDimension(vtkType);
805   if (cellDim != 3)
806     return 0; // TODO voisins des faces ou edges
807   int cellId = this->_cellIdToDownId[vtkId];
808
809   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
810   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
811   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
812
813   // --- iteration on faces of the 3D cell.
814
815   int nb = 0;
816   for (int i = 0; i < nbCells; i++)
817     {
818       int downId = downCells[i];
819       int cellType = downTyp[i];
820       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
821       const int *upCells = _downArray[cellType]->getUpCells(downId);
822       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
823
824       // --- max 2 upCells, one is this cell, the other is a neighbor
825
826       for (int j = 0; j < nbUp; j++)
827         {
828           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
829             continue;
830           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
831           neighborsVtkIds[nb] = vtkNeighbor;
832           downIds[nb] = downId;
833           downTypes[nb] = cellType;
834           nb++;
835         }
836       if (nb >= NBMAXNEIGHBORS)
837         assert(0);
838     }
839   return nb;
840 }
841
842 /*! get the node id's of a cell.
843  * The cell is defined by it's downward connectivity id and type.
844  * @param nodeSet set of of vtk node id's to fill.
845  * @param downId downward connectivity id of the cell.
846  * @param downType type of cell.
847  */
848 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
849 {
850   _downArray[downType]->getNodeIds(downId, nodeSet);
851 }
852
853 /*! change some nodes in cell without modifying type or internal connectivity.
854  * Nodes inverse connectivity is maintained up to date.
855  * @param vtkVolId vtk id of the cell
856  * @param localClonedNodeIds map old node id to new node id.
857  */
858 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
859 {
860   vtkIdType npts = 0;
861   vtkIdType *pts; // will refer to the point id's of the face
862   this->GetCellPoints(vtkVolId, npts, pts);
863   for (int i = 0; i < npts; i++)
864     {
865       if (localClonedNodeIds.count(pts[i]))
866         {
867           vtkIdType oldpt = pts[i];
868           pts[i] = localClonedNodeIds[oldpt];
869           //MESSAGE(oldpt << " --> " << pts[i]);
870           //this->RemoveReferenceToCell(oldpt, vtkVolId);
871           //this->AddReferenceToCell(pts[i], vtkVolId);
872         }
873     }
874 }
875
876 /*! Create a volume (prism or hexahedron) by duplication of a face.
877  * the nodes of the new face are already created.
878  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
879  * @param localClonedNodeIds map old node id to new node id.
880  * @return vtk id of the new volume.
881  */
882 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<int>& orderedNodes)
883 {
884   int vtkType = this->GetCellType(vtkVolId);
885   int cellDim = SMDS_Downward::getCellDimension(vtkType);
886   if (cellDim != 3)
887     return 0;
888   SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
889   int downVolId = this->_cellIdToDownId[vtkVolId];
890   downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
891   return orderedNodes.size();
892 }
893
894 void SMDS_UnstructuredGrid::BuildLinks()
895 {
896   // Remove the old links if they are already built
897   if (this->Links)
898     {
899     this->Links->UnRegister(this);
900     }
901
902   this->Links = SMDS_CellLinks::New();
903   this->Links->Allocate(this->GetNumberOfPoints());
904   this->Links->Register(this);
905   this->Links->BuildLinks(this, this->Connectivity);
906   this->Links->Delete();
907 }