Salome HOME
a7837d4934857ced24809dfa38aaab838c339938
[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   MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
135   int alreadyCopied = 0;
136
137   // --- if newNodeSize, create a new compacted vtkPoints
138
139   vtkPoints *newPoints = vtkPoints::New();
140   newPoints->SetDataType(VTK_DOUBLE);
141   newPoints->SetNumberOfPoints(newNodeSize);
142   if (newNodeSize)
143     {
144       MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
145       // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
146       // using double type for storing coordinates of nodes instead float.
147       int oldNodeSize = idNodesOldToNew.size();
148
149       int i = 0;
150       while ( i < oldNodeSize )
151       {
152         // skip a hole if any
153         while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
154           ++i;
155         int startBloc = i;
156         // look for a block end
157         while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
158           ++i;
159         int endBloc = i;
160         copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
161       }
162       newPoints->Squeeze();
163     }
164
165   // --- create new compacted Connectivity, Locations and Types
166
167   int oldCellSize = this->Types->GetNumberOfTuples();
168
169   vtkCellArray *newConnectivity = vtkCellArray::New();
170   newConnectivity->Initialize();
171   int oldCellDataSize = this->Connectivity->GetData()->GetSize();
172   newConnectivity->Allocate(oldCellDataSize);
173   MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
174
175   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
176   newTypes->Initialize();
177   newTypes->SetNumberOfValues(newCellSize);
178
179   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
180   newLocations->Initialize();
181   newLocations->SetNumberOfValues(newCellSize);
182
183   // TODO some polyhedron may be huge (only in some tests)
184   vtkIdType tmpid[NBMAXNODESINCELL];
185   vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
186
187   alreadyCopied = 0;
188   int i = 0;
189   while ( i < oldCellSize )
190   {
191     // skip a hole if any
192     while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
193       ++i;
194     int startBloc = i;
195     // look for a block end
196     while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
197       ++i;
198     int endBloc = i;
199     if ( endBloc > startBloc )
200       copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
201                alreadyCopied, startBloc, endBloc);
202   }
203
204   newConnectivity->Squeeze();
205
206   if (1/*newNodeSize*/)
207     {
208       MESSAGE("------- newNodeSize, setPoints");
209       this->SetPoints(newPoints);
210       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
211     }
212
213   if (this->FaceLocations)
214     {
215       vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
216       newFaceLocations->Initialize();
217       newFaceLocations->Allocate(newTypes->GetSize());
218       vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
219       newFaces->Initialize();
220       newFaces->Allocate(this->Faces->GetSize());
221       for (int i = 0; i < oldCellSize; i++)
222         {
223           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
224             continue;
225           int newCellId = idCellsOldToNew[i];
226           if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
227             {
228                newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
229                int oldFaceLoc = this->FaceLocations->GetValue(i);
230                int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
231                newFaces->InsertNextValue(nCellFaces);
232                for (int n=0; n<nCellFaces; n++)
233                  {
234                    int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
235                    newFaces->InsertNextValue(nptsInFace);
236                    for (int k=0; k<nptsInFace; k++)
237                      {
238                        int oldpt = this->Faces->GetValue(oldFaceLoc++);
239                        newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
240                      }
241                  }
242             }
243           else
244             {
245                newFaceLocations->InsertNextValue(-1);
246             }
247         }
248       newFaceLocations->Squeeze();
249       newFaces->Squeeze();
250       newFaceLocations->Register(this);
251       newFaces->Register(this);
252       this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
253       newFaceLocations->Delete();
254       newFaces->Delete();
255     }
256   else
257     this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
258
259   newPoints->Delete();
260   newTypes->Delete();
261   newLocations->Delete();
262   newConnectivity->Delete();
263   this->BuildLinks();
264 }
265
266 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
267                                       int start, int end)
268 {
269   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
270   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
271   void *source = this->Points->GetVoidPointer(3 * start);
272   int nbPoints = end - start;
273   if (nbPoints > 0)
274     {
275       memcpy(target, source, 3 * sizeof(double) * nbPoints);
276       for (int j = start; j < end; j++)
277         idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
278     }
279 }
280
281 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
282                                      std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
283                                      vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
284                                      int start, int end)
285 {
286   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
287   for (int j = start; j < end; j++)
288     {
289       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
290       idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
291       vtkIdType oldLoc = this->Locations->GetValue(j);
292       vtkIdType nbpts;
293       vtkIdType *oldPtsCell = 0;
294       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
295       assert(nbpts < NBMAXNODESINCELL);
296       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
297       for (int l = 0; l < nbpts; l++)
298         {
299           int oldval = oldPtsCell[l];
300           pointsCell[l] = idNodesOldToNew[oldval];
301           //MESSAGE("   " << oldval << " " << pointsCell[l]);
302         }
303       /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
304       int newLoc = newConnectivity->GetInsertLocation(nbpts);
305       //MESSAGE(newcnt << " " << newLoc);
306       newLocations->SetValue(alreadyCopied, newLoc);
307       alreadyCopied++;
308     }
309 }
310
311 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
312 {
313   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
314   return _cellIdToDownId[vtkCellId];
315 }
316
317 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
318 {
319   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
320   _cellIdToDownId[vtkCellId] = downId;
321 }
322
323 /*! Build downward connectivity: to do only when needed because heavy memory load.
324  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
325  *
326  */
327 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
328 {
329   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
330   // TODO calcul partiel sans edges
331
332   // --- erase previous data if any
333
334   for (int i = 0; i < _downArray.size(); i++)
335     {
336       if (_downArray[i])
337         delete _downArray[i];
338       _downArray[i] = 0;
339     }
340   _cellIdToDownId.clear();
341
342   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
343
344   _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
345
346   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
347   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
348   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
349   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
350   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
351   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
352   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
353   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
354   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
355   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
356   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
357   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
358   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
359   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
360
361   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
362
363   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
364
365   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
366   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
367   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
368   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
369   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
370   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
371   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
372   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
373
374   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
375   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
376   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
377   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
378   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
379   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
380
381   int GuessSize[VTK_QUADRATIC_TETRA];
382   GuessSize[VTK_LINE] = nbLineGuess;
383   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
384   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
385   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
386   GuessSize[VTK_QUAD] = nbLinQuadGuess;
387   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
388   GuessSize[VTK_TETRA] = nbLinTetra;
389   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
390   GuessSize[VTK_PYRAMID] = nbLinPyra;
391   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
392   GuessSize[VTK_WEDGE] = nbLinPrism;
393   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
394   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
395   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
396
397   _downArray[VTK_LINE]->allocate(nbLineGuess);
398   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
399   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
400   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
401   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
402   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
403   _downArray[VTK_TETRA]->allocate(nbLinTetra);
404   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
405   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
406   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
407   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
408   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
409   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
410   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
411
412   // --- iteration on vtkUnstructuredGrid cells, only faces
413   //     for each vtk face:
414   //       create a downward face entry with its downward id.
415   //       compute vtk volumes, create downward volumes entry.
416   //       mark face in downward volumes
417   //       mark volumes in downward face
418
419   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
420   int cellSize = this->Types->GetNumberOfTuples();
421   _cellIdToDownId.resize(cellSize, -1);
422
423   for (int i = 0; i < cellSize; i++)
424     {
425       int vtkFaceType = this->GetCellType(i);
426       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
427         {
428           int vtkFaceId = i;
429           //ASSERT(_downArray[vtkFaceType]);
430           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
431           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
432           downFace->setTempNodes(connFaceId, vtkFaceId);
433           int vols[2] = { -1, -1 };
434           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
435           //MESSAGE("nbVolumes="<< nbVolumes);
436           for (int ivol = 0; ivol < nbVolumes; ivol++)
437             {
438               int vtkVolId = vols[ivol];
439               int vtkVolType = this->GetCellType(vtkVolId);
440               //ASSERT(_downArray[vtkVolType]);
441               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
442               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
443               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
444               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
445             }
446         }
447     }
448
449   // --- iteration on vtkUnstructuredGrid cells, only volumes
450   //     for each vtk volume:
451   //       create downward volumes entry if not already done
452   //       build a temporary list of faces described with their nodes
453   //       for each face
454   //         compute the vtk volumes containing this face
455   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
456   //         if not, create a downward face entry (resizing of structure required)
457   //         (the downward faces store a temporary list of nodes to ease the comparison)
458   //         create downward volumes entry if not already done
459   //         mark volumes in downward face
460   //         mark face in downward volumes
461
462   CHRONOSTOP(20);
463   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
464
465   for (int i = 0; i < cellSize; i++)
466     {
467       int vtkType = this->GetCellType(i);
468       if (SMDS_Downward::getCellDimension(vtkType) == 3)
469         {
470           //CHRONO(31);
471           int vtkVolId = i;
472           // MESSAGE("vtk volume " << vtkVolId);
473           //ASSERT(_downArray[vtkType]);
474           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
475
476           // --- find all the faces of the volume, describe the faces by their nodes
477
478           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
479           ListElemByNodesType facesWithNodes;
480           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
481           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
482           //CHRONOSTOP(31);
483           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
484             {
485               // --- find the volumes containing the face
486
487               //CHRONO(32);
488               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
489               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
490               int vols[2] = { -1, -1 };
491               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
492               int lg = facesWithNodes.elems[iface].nbNodes;
493               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
494               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
495
496               // --- check if face is registered in the volumes
497               //CHRONOSTOP(32);
498
499               //CHRONO(33);
500               int connFaceId = -1;
501               for (int ivol = 0; ivol < nbVolumes; ivol++)
502                 {
503                   int vtkVolId2 = vols[ivol];
504                   int vtkVolType = this->GetCellType(vtkVolId2);
505                   //ASSERT(_downArray[vtkVolType]);
506                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
507                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
508                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
509                   if (connFaceId >= 0)
510                     break; // --- face already created
511                 }//CHRONOSTOP(33);
512
513               // --- if face is not registered in the volumes, create face
514
515               //CHRONO(34);
516               if (connFaceId < 0)
517                 {
518                   connFaceId = _downArray[vtkFaceType]->addCell();
519                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
520                 }//CHRONOSTOP(34);
521
522               // --- mark volumes in downward face and mark face in downward volumes
523
524               //CHRONO(35);
525               for (int ivol = 0; ivol < nbVolumes; ivol++)
526                 {
527                   int vtkVolId2 = vols[ivol];
528                   int vtkVolType = this->GetCellType(vtkVolId2);
529                   //ASSERT(_downArray[vtkVolType]);
530                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
531                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
532                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
533                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
534                 }//CHRONOSTOP(35);
535             }
536         }
537     }
538
539   // --- iteration on vtkUnstructuredGrid cells, only edges
540   //     for each vtk edge:
541   //       create downward edge entry
542   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
543   //       find downward faces
544   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
545   //       mark edge in downward faces
546   //       mark faces in downward edge
547
548   CHRONOSTOP(21);
549   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
550
551   for (int i = 0; i < cellSize; i++)
552     {
553       int vtkEdgeType = this->GetCellType(i);
554       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
555         {
556           int vtkEdgeId = i;
557           //ASSERT(_downArray[vtkEdgeType]);
558           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
559           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
560           downEdge->setNodes(connEdgeId, vtkEdgeId);
561           vector<int> vtkIds;
562           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
563           int downFaces[1000];
564           unsigned char downTypes[1000];
565           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
566           for (int n = 0; n < nbDownFaces; n++)
567             {
568               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
569               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
570             }
571         }
572     }
573
574   // --- iteration on downward faces (they are all listed now)
575   //     for each downward face:
576   //       build a temporary list of edges with their ordered list of nodes
577   //       for each edge:
578   //         find all the vtk cells containing this edge
579   //         then identify all the downward faces containing the edge, from the vtk cells
580   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
581   //         if not, create a downward edge entry with the node id's
582   //         mark edge in downward faces
583   //         mark downward faces in edge (size of list unknown, to be allocated)
584
585   CHRONOSTOP(22);CHRONO(23);
586
587   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
588     {
589       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
590         continue;
591
592       // --- find all the edges of the face, describe the edges by their nodes
593
594       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
595       int maxId = downFace->getMaxId();
596       for (int faceId = 0; faceId < maxId; faceId++)
597         {
598           //CHRONO(40);
599           ListElemByNodesType edgesWithNodes;
600           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
601           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
602
603           //CHRONOSTOP(40);
604           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
605             {
606
607               // --- check if the edge is already registered by exploration of the faces
608
609               //CHRONO(41);
610               vector<int> vtkIds;
611               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
612               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
613               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
614               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
615               //CHRONOSTOP(41);CHRONO(42);
616               int downFaces[1000];
617               unsigned char downTypes[1000];
618               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
619               //CHRONOSTOP(42);
620
621               //CHRONO(43);
622               int connEdgeId = -1;
623               for (int idf = 0; idf < nbDownFaces; idf++)
624                 {
625                   int faceId2 = downFaces[idf];
626                   int faceType = downTypes[idf];
627                   //ASSERT(_downArray[faceType]);
628                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
629                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
630                   if (connEdgeId >= 0)
631                     break; // --- edge already created
632                 }//CHRONOSTOP(43);
633
634               // --- if edge is not registered in the faces, create edge
635
636               if (connEdgeId < 0)
637                 {
638                   //CHRONO(44);
639                   connEdgeId = _downArray[vtkEdgeType]->addCell();
640                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
641                   //CHRONOSTOP(44);
642                 }
643
644               // --- mark faces in downward edge and mark edge in downward faces
645
646               //CHRONO(45);
647               for (int idf = 0; idf < nbDownFaces; idf++)
648                 {
649                   int faceId2 = downFaces[idf];
650                   int faceType = downTypes[idf];
651                   //ASSERT(_downArray[faceType]);
652                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
653                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
654                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
655                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
656                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
657                 }//CHRONOSTOP(45);
658             }
659         }
660     }
661
662   CHRONOSTOP(23);CHRONO(24);
663
664   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
665   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
666
667   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
668     {
669       if (SMDS_Downward *down = _downArray[vtkType])
670         {
671           down->compactStorage();
672         }
673     }
674
675   // --- Statistics
676
677   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
678     {
679       if (SMDS_Downward *down = _downArray[vtkType])
680         {
681           if (down->getMaxId())
682             {
683               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
684                   << GuessSize[vtkType] << " real: " << down->getMaxId());
685             }
686         }
687     }CHRONOSTOP(24);CHRONOSTOP(2);
688   counters::stats();
689 }
690
691 /*! Get the neighbors of a cell.
692  * Only the neighbors having the dimension of the cell are taken into account
693  * (neighbors of a volume are the volumes sharing a face with this volume,
694  *  neighbors of a face are the faces sharing an edge with this face...).
695  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
696  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
697  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
698  * @param vtkId the vtk id of the cell
699  * @return number of neighbors
700  */
701 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
702 {
703   int vtkType = this->GetCellType(vtkId);
704   int cellDim = SMDS_Downward::getCellDimension(vtkType);
705   if (cellDim <2)
706     return 0; // TODO voisins des edges = edges connectees
707   int cellId = this->_cellIdToDownId[vtkId];
708
709   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
710   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
711   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
712
713   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
714
715   int nb = 0;
716   for (int i = 0; i < nbCells; i++)
717     {
718       int downId = downCells[i];
719       int cellType = downTyp[i];
720       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
721       const int *upCells = _downArray[cellType]->getUpCells(downId);
722       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
723
724       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
725       //    for a face, number of neighbors (connected faces) not known
726
727       for (int j = 0; j < nbUp; j++)
728         {
729           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
730             continue;
731           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
732           neighborsVtkIds[nb] = vtkNeighbor;
733           downIds[nb] = downId;
734           downTypes[nb] = cellType;
735           nb++;
736         }
737       if (nb >= NBMAXNEIGHBORS)
738         assert(0);
739     }
740   return nb;
741 }
742
743 /*! get the node id's of a cell.
744  * The cell is defined by it's downward connectivity id and type.
745  * @param nodeSet set of of vtk node id's to fill.
746  * @param downId downward connectivity id of the cell.
747  * @param downType type of cell.
748  */
749 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
750 {
751   _downArray[downType]->getNodeIds(downId, nodeSet);
752 }
753
754 /*! change some nodes in cell without modifying type or internal connectivity.
755  * Nodes inverse connectivity is maintained up to date.
756  * @param vtkVolId vtk id of the cell
757  * @param localClonedNodeIds map old node id to new node id.
758  */
759 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
760 {
761   vtkIdType npts = 0;
762   vtkIdType *pts; // will refer to the point id's of the face
763   this->GetCellPoints(vtkVolId, npts, pts);
764   for (int i = 0; i < npts; i++)
765     {
766       if (localClonedNodeIds.count(pts[i]))
767         {
768           vtkIdType oldpt = pts[i];
769           pts[i] = localClonedNodeIds[oldpt];
770           //MESSAGE(oldpt << " --> " << pts[i]);
771           //this->RemoveReferenceToCell(oldpt, vtkVolId);
772           //this->AddReferenceToCell(pts[i], vtkVolId);
773         }
774     }
775 }
776
777 /*! reorder the nodes of a face
778  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
779  * @param orderedNodes list of nodes to reorder (in out)
780  * @return size of the list
781  */
782 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes)
783 {
784   int vtkType = this->GetCellType(vtkVolId);
785   int cellDim = SMDS_Downward::getCellDimension(vtkType);
786   if (cellDim != 3)
787     return 0;
788   SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
789   int downVolId = this->_cellIdToDownId[vtkVolId];
790   downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
791   return orderedNodes.size();
792 }
793
794 void SMDS_UnstructuredGrid::BuildLinks()
795 {
796   // Remove the old links if they are already built
797   if (this->Links)
798     {
799     this->Links->UnRegister(this);
800     }
801
802   this->Links = SMDS_CellLinks::New();
803   this->Links->Allocate(this->GetNumberOfPoints());
804   this->Links->Register(this);
805   this->Links->BuildLinks(this, this->Connectivity);
806   this->Links->Delete();
807 }