Salome HOME
merge V5_1_4
[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_UnstructuredGrid* SMDS_UnstructuredGrid::New()
19 {
20   MESSAGE("SMDS_UnstructuredGrid::New");
21   return new SMDS_UnstructuredGrid();
22 }
23
24 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
25   vtkUnstructuredGrid()
26 {
27   _cellIdToDownId.clear();
28   _downTypes.clear();
29   _downArray.clear();
30   _mesh = 0;
31   _counters = new counters(100);
32 }
33
34 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
35 {
36 }
37
38 unsigned long SMDS_UnstructuredGrid::GetMTime()
39 {
40   unsigned long mtime = vtkUnstructuredGrid::GetMTime();
41   MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
42   return mtime;
43 }
44
45 void SMDS_UnstructuredGrid::Update()
46 {
47   MESSAGE("SMDS_UnstructuredGrid::Update");
48   return vtkUnstructuredGrid::Update();
49 }
50
51 void SMDS_UnstructuredGrid::UpdateInformation()
52 {
53   MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
54   return vtkUnstructuredGrid::UpdateInformation();
55 }
56
57 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
58 {
59   _mesh = mesh;
60 }
61
62 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
63                                         std::vector<int>& idCellsOldToNew, int newCellSize)
64 {
65   MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
66   int startHole = 0;
67   int endHole = 0;
68   int startBloc = 0;
69   int endBloc = 0;
70   int alreadyCopied = 0;
71   int holes = 0;
72
73   typedef enum
74   {
75     lookHoleStart, lookHoleEnd, lookBlocEnd
76   } enumState;
77   enumState compactState = lookHoleStart;
78
79   //    if (this->Links)
80   //    {
81   //            this->Links->UnRegister(this);
82   //            this->Links = 0;
83   //    }
84
85   // --- if newNodeSize, create a new compacted vtkPoints
86
87   vtkPoints *newPoints = 0;
88   if (newNodeSize)
89     {
90       MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
91       newPoints = vtkPoints::New();
92       newPoints->Initialize();
93       newPoints->Allocate(newNodeSize);
94       newPoints->SetNumberOfPoints(newNodeSize);
95       int oldNodeSize = idNodesOldToNew.size();
96
97       for (int i = 0; i < oldNodeSize; i++)
98         {
99           switch (compactState)
100           {
101             case lookHoleStart:
102               if (idNodesOldToNew[i] < 0)
103                 {
104                   MESSAGE("-------------- newNodeSize, startHole " << i << " " << oldNodeSize);
105                   startHole = i;
106                   if (!alreadyCopied) // copy the first bloc
107                     {
108                       MESSAGE("--------- copy first nodes before hole " << i << " " << oldNodeSize);
109                       copyNodes(newPoints, idNodesOldToNew, alreadyCopied, 0, startHole);
110                     }
111                   compactState = lookHoleEnd;
112                 }
113               break;
114             case lookHoleEnd:
115               if (idNodesOldToNew[i] >= 0)
116                 {
117                   MESSAGE("-------------- newNodeSize, endHole " << i << " " << oldNodeSize);
118                   endHole = i;
119                   startBloc = i;
120                   compactState = lookBlocEnd;
121                 }
122               break;
123             case lookBlocEnd:
124               if (idNodesOldToNew[i] < 0)
125                 endBloc = i; // see nbPoints below
126               else if (i == (oldNodeSize - 1))
127                 endBloc = i + 1;
128               if (endBloc)
129                 {
130                   MESSAGE("-------------- newNodeSize, endbloc " << endBloc << " " << oldNodeSize);
131                   copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
132                   compactState = lookHoleStart;
133                   startHole = i;
134                   endHole = 0;
135                   startBloc = 0;
136                   endBloc = 0;
137                 }
138               break;
139           }
140         }
141       if (!alreadyCopied) // no hole, but shorter, no need to modify idNodesOldToNew
142         {
143           MESSAGE("------------- newNodeSize, shorter " << oldNodeSize);
144           copyNodes(newPoints, idNodesOldToNew, alreadyCopied, 0, newNodeSize);
145         }
146       newPoints->Squeeze();
147     }
148
149   // --- create new compacted Connectivity, Locations and Types
150
151   int oldCellSize = this->Types->GetNumberOfTuples();
152
153   vtkCellArray *newConnectivity = vtkCellArray::New();
154   newConnectivity->Initialize();
155   int oldCellDataSize = this->Connectivity->GetData()->GetSize();
156   newConnectivity->Allocate(oldCellDataSize);
157   MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
158
159   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
160   newTypes->Initialize();
161   //newTypes->Allocate(oldCellSize);
162   newTypes->SetNumberOfValues(newCellSize);
163
164   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
165   newLocations->Initialize();
166   //newLocations->Allocate(oldCellSize);
167   newLocations->SetNumberOfValues(newCellSize);
168
169   startHole = 0;
170   endHole = 0;
171   startBloc = 0;
172   endBloc = 0;
173   alreadyCopied = 0;
174   holes = 0;
175   compactState = lookHoleStart;
176
177   vtkIdType tmpid[50];
178   vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
179
180   for (int i = 0; i < oldCellSize; i++)
181     {
182       switch (compactState)
183       {
184         case lookHoleStart:
185           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
186             {
187               MESSAGE(" -------- newCellSize, startHole " << i << " " << oldCellSize);
188               startHole = i;
189               compactState = lookHoleEnd;
190               if (!alreadyCopied) // copy the first bloc
191                 {
192                   MESSAGE("--------- copy first bloc before hole " << i << " " << oldCellSize);
193                   copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
194                            alreadyCopied, 0, startHole);
195                 }
196             }
197           break;
198         case lookHoleEnd:
199           if (this->Types->GetValue(i) != VTK_EMPTY_CELL)
200             {
201               MESSAGE(" -------- newCellSize, EndHole " << i << " " << oldCellSize);
202               endHole = i;
203               startBloc = i;
204               compactState = lookBlocEnd;
205               holes += endHole - startHole;
206               //alreadyCopied = startBloc -holes;
207             }
208           break;
209         case lookBlocEnd:
210           endBloc = 0;
211           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
212             endBloc = i;
213           else if (i == (oldCellSize - 1))
214             endBloc = i + 1;
215           if (endBloc)
216             {
217               MESSAGE(" -------- newCellSize, endBloc " << endBloc << " " << oldCellSize);
218               copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
219                        alreadyCopied, startBloc, endBloc);
220               compactState = lookHoleStart;
221             }
222           break;
223       }
224     }
225   if (!alreadyCopied) // no hole, but shorter
226     {
227       MESSAGE(" -------- newCellSize, shorter " << oldCellSize);
228       copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell, alreadyCopied, 0,
229                oldCellSize);
230     }
231
232   newConnectivity->Squeeze();
233   //newTypes->Squeeze();
234   //newLocations->Squeeze();
235
236   if (newNodeSize)
237     {
238       MESSAGE("------- newNodeSize, setPoints");
239       this->SetPoints(newPoints);
240       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
241     }
242   this->SetCells(newTypes, newLocations, newConnectivity);
243   this->BuildLinks();
244 }
245
246 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
247                                       int start, int end)
248 {
249   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
250   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
251   void *source = this->Points->GetVoidPointer(3 * start);
252   int nbPoints = end - start;
253   if (nbPoints > 0)
254     {
255       memcpy(target, source, 3 * sizeof(float) * nbPoints);
256       for (int j = start; j < end; j++)
257         idNodesOldToNew[j] = alreadyCopied++;
258     }
259 }
260
261 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
262                                      std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
263                                      vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
264                                      int start, int end)
265 {
266   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
267   for (int j = start; j < end; j++)
268     {
269       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
270       idCellsOldToNew[j] = alreadyCopied;
271       vtkIdType oldLoc = this->Locations->GetValue(j);
272       vtkIdType nbpts;
273       vtkIdType *oldPtsCell = 0;
274       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
275       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
276       for (int l = 0; l < nbpts; l++)
277         {
278           int oldval = oldPtsCell[l];
279           pointsCell[l] = idNodesOldToNew[oldval];
280           //MESSAGE("   " << oldval << " " << pointsCell[l]);
281         }
282       int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
283       int newLoc = newConnectivity->GetInsertLocation(nbpts);
284       //MESSAGE(newcnt << " " << newLoc);
285       newLocations->SetValue(alreadyCopied, newLoc);
286       alreadyCopied++;
287     }
288 }
289
290 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
291 {
292   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
293   return _cellIdToDownId[vtkCellId];
294 }
295
296 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
297 {
298   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
299   _cellIdToDownId[vtkCellId] = downId;
300 }
301
302 /*! Build downward connectivity: to do only when needed because heavy memory load.
303  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
304  *
305  */
306 void SMDS_UnstructuredGrid::BuildDownwardConnectivity()
307 {
308   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
309
310   // --- erase previous data if any
311
312   for (int i = 0; i < _downArray.size(); i++)
313     {
314       if (_downArray[i])
315         delete _downArray[i];
316       _downArray[i] = 0;
317     }
318   _cellIdToDownId.clear();
319
320   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
321
322   _downArray.resize(VTK_QUADRATIC_PYRAMID + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
323
324   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
325   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
326   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
327   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
328   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
329   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
330   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
331   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
332   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
333   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
334   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
335   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
336   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
337   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
338
339   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
340
341   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
342
343   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
344   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
345   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
346   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
347   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
348   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
349   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
350   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
351
352   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
353   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
354   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
355   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
356   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
357   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
358
359   int GuessSize[VTK_QUADRATIC_TETRA];
360   GuessSize[VTK_LINE] = nbLineGuess;
361   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
362   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
363   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
364   GuessSize[VTK_QUAD] = nbLinQuadGuess;
365   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
366   GuessSize[VTK_TETRA] = nbLinTetra;
367   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
368   GuessSize[VTK_PYRAMID] = nbLinPyra;
369   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
370   GuessSize[VTK_WEDGE] = nbLinPrism;
371   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
372   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
373   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
374
375   _downArray[VTK_LINE]->allocate(nbLineGuess);
376   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
377   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
378   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
379   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
380   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
381   _downArray[VTK_TETRA]->allocate(nbLinTetra);
382   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
383   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
384   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
385   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
386   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
387   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
388   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
389
390   // --- iteration on vtkUnstructuredGrid cells, only faces
391   //     for each vtk face:
392   //       create a downward face entry with its downward id.
393   //       compute vtk volumes, create downward volumes entry.
394   //       mark face in downward volumes
395   //       mark volumes in downward face
396
397   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
398   int cellSize = this->Types->GetNumberOfTuples();
399   _cellIdToDownId.resize(cellSize, -1);
400
401   for (int i = 0; i < cellSize; i++)
402     {
403       int vtkFaceType = this->GetCellType(i);
404       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
405         {
406           int vtkFaceId = i;
407           //ASSERT(_downArray[vtkFaceType]);
408           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
409           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
410           downFace->setTempNodes(connFaceId, vtkFaceId);
411           int vols[2] = { -1, -1 };
412           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
413           //MESSAGE("nbVolumes="<< nbVolumes);
414           for (int ivol = 0; ivol < nbVolumes; ivol++)
415             {
416               int vtkVolId = vols[ivol];
417               int vtkVolType = this->GetCellType(vtkVolId);
418               //ASSERT(_downArray[vtkVolType]);
419               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
420               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
421               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
422               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
423             }
424         }
425     }
426
427   // --- iteration on vtkUnstructuredGrid cells, only volumes
428   //     for each vtk volume:
429   //       create downward volumes entry if not already done
430   //       build a temporary list of faces described with their nodes
431   //       for each face
432   //         compute the vtk volumes containing this face
433   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
434   //         if not, create a downward face entry (resizing of structure required)
435   //         (the downward faces store a temporary list of nodes to ease the comparison)
436   //         create downward volumes entry if not already done
437   //         mark volumes in downward face
438   //         mark face in downward volumes
439
440   CHRONOSTOP(20);
441   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
442
443   for (int i = 0; i < cellSize; i++)
444     {
445       int vtkType = this->GetCellType(i);
446       if (SMDS_Downward::getCellDimension(vtkType) == 3)
447         {
448           //CHRONO(31);
449           int vtkVolId = i;
450           // MESSAGE("vtk volume " << vtkVolId);
451           //ASSERT(_downArray[vtkType]);
452           int connVolId = _downArray[vtkType]->addCell(vtkVolId);
453
454           // --- find all the faces of the volume, describe the faces by their nodes
455
456           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
457           ListElemByNodesType facesWithNodes;
458           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
459           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
460           //CHRONOSTOP(31);
461           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
462             {
463               // --- find the volumes containing the face
464
465               //CHRONO(32);
466               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
467               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
468               int vols[2] = { -1, -1 };
469               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
470               int lg = facesWithNodes.elems[iface].nbNodes;
471               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
472               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
473
474               // --- check if face is registered in the volumes
475               //CHRONOSTOP(32);
476
477               //CHRONO(33);
478               int connFaceId = -1;
479               for (int ivol = 0; ivol < nbVolumes; ivol++)
480                 {
481                   int vtkVolId2 = vols[ivol];
482                   int vtkVolType = this->GetCellType(vtkVolId2);
483                   //ASSERT(_downArray[vtkVolType]);
484                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
485                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
486                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
487                   if (connFaceId >= 0)
488                     break; // --- face already created
489                 }//CHRONOSTOP(33);
490
491               // --- if face is not registered in the volumes, create face
492
493               //CHRONO(34);
494               if (connFaceId < 0)
495                 {
496                   connFaceId = _downArray[vtkFaceType]->addCell();
497                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
498                 }//CHRONOSTOP(34);
499
500               // --- mark volumes in downward face and mark face in downward volumes
501
502               //CHRONO(35);
503               for (int ivol = 0; ivol < nbVolumes; ivol++)
504                 {
505                   int vtkVolId2 = vols[ivol];
506                   int vtkVolType = this->GetCellType(vtkVolId2);
507                   //ASSERT(_downArray[vtkVolType]);
508                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
509                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
510                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
511                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
512                 }//CHRONOSTOP(35);
513             }
514         }
515     }
516
517   // --- iteration on vtkUnstructuredGrid cells, only edges
518   //     for each vtk edge:
519   //       create downward edge entry
520   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
521   //       find downward faces
522   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
523   //       mark edge in downward faces
524   //       mark faces in downward edge
525
526   CHRONOSTOP(21);
527   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
528
529   for (int i = 0; i < cellSize; i++)
530     {
531       int vtkEdgeType = this->GetCellType(i);
532       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
533         {
534           int vtkEdgeId = i;
535           //ASSERT(_downArray[vtkEdgeType]);
536           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
537           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
538           downEdge->setNodes(connEdgeId, vtkEdgeId);
539           vector<int> vtkIds;
540           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
541           int downFaces[1000];
542           unsigned char downTypes[1000];
543           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
544           for (int n = 0; n < nbDownFaces; n++)
545             {
546               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
547               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
548             }
549         }
550     }
551
552   // --- iteration on downward faces (they are all listed now)
553   //     for each downward face:
554   //       build a temporary list of edges with their ordered list of nodes
555   //       for each edge:
556   //         find all the vtk cells containing this edge
557   //         then identify all the downward faces containing the edge, from the vtk cells
558   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
559   //         if not, create a downward edge entry with the node id's
560   //         mark edge in downward faces
561   //         mark downward faces in edge (size of list unknown, to be allocated)
562
563   CHRONOSTOP(22);CHRONO(23);
564
565   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
566     {
567       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
568         continue;
569
570       // --- find all the edges of the face, describe the edges by their nodes
571
572       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
573       int maxId = downFace->getMaxId();
574       for (int faceId = 0; faceId < maxId; faceId++)
575         {
576           //CHRONO(40);
577           ListElemByNodesType edgesWithNodes;
578           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
579           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
580
581           //CHRONOSTOP(40);
582           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
583             {
584
585               // --- check if the edge is already registered by exploration of the faces
586
587               //CHRONO(41);
588               vector<int> vtkIds;
589               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
590               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
591               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
592               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
593               //CHRONOSTOP(41);CHRONO(42);
594               int downFaces[1000];
595               unsigned char downTypes[1000];
596               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
597               //CHRONOSTOP(42);
598
599               //CHRONO(43);
600               int connEdgeId = -1;
601               for (int idf = 0; idf < nbDownFaces; idf++)
602                 {
603                   int faceId2 = downFaces[idf];
604                   int faceType = downTypes[idf];
605                   //ASSERT(_downArray[faceType]);
606                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
607                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
608                   if (connEdgeId >= 0)
609                     break; // --- edge already created
610                 }//CHRONOSTOP(43);
611
612               // --- if edge is not registered in the faces, create edge
613
614               if (connEdgeId < 0)
615                 {
616                   //CHRONO(44);
617                   connEdgeId = _downArray[vtkEdgeType]->addCell();
618                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
619                   //CHRONOSTOP(44);
620                 }
621
622               // --- mark faces in downward edge and mark edge in downward faces
623
624               //CHRONO(45);
625               for (int idf = 0; idf < nbDownFaces; idf++)
626                 {
627                   int faceId2 = downFaces[idf];
628                   int faceType = downTypes[idf];
629                   //ASSERT(_downArray[faceType]);
630                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
631                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
632                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
633                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
634                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
635                 }//CHRONOSTOP(45);
636             }
637         }
638     }
639
640   CHRONOSTOP(23);CHRONO(24);
641
642   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
643   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
644
645   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
646     {
647       if (SMDS_Downward *down = _downArray[vtkType])
648         {
649           down->compactStorage();
650         }
651     }
652
653   // --- Statistics
654
655   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
656     {
657       if (SMDS_Downward *down = _downArray[vtkType])
658         {
659           if (down->getMaxId())
660             {
661               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
662                   << GuessSize[vtkType] << " real: " << down->getMaxId());
663             }
664         }
665     }CHRONOSTOP(24);CHRONOSTOP(2);
666   _counters->stats();
667 }