Salome HOME
PR: correct some warnings
[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 //#ifdef VTK_HAVE_POLYHEDRON
313   // TODO compact faces for Polyhedrons
314   // refaire completement faces et faceLocation
315   // pour chaque cell, recup oldCellId, oldFacesId, recopie dans newFaces de la faceStream
316   // en changeant les numeros de noeuds
317 //  vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
318 //  newFaceLocations->Initialize();
319 //  vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
320 //  newFaces->Initialize();
321 //  newFaceLocations->DeepCopy(this->FaceLocations);
322 //  newFaces->DeepCopy(this->Faces);
323 //  this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
324 //  newFaceLocations->Delete();
325 //  newFaces->Delete();
326   if (this->FaceLocations) this->FaceLocations->Register(this);
327   if (this->Faces) this->Faces->Register(this);
328   this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
329 //#else
330 //  this->SetCells(newTypes, newLocations, newConnectivity);
331 //#endif
332   newTypes->Delete();
333   newLocations->Delete();
334   newConnectivity->Delete();
335   this->BuildLinks();
336 }
337
338 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
339                                       int start, int end)
340 {
341   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
342   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
343   void *source = this->Points->GetVoidPointer(3 * start);
344   int nbPoints = end - start;
345   if (nbPoints > 0)
346     {
347       memcpy(target, source, 3 * sizeof(float) * nbPoints);
348       for (int j = start; j < end; j++)
349         idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
350         //idNodesOldToNew[alreadyCopied++] = idNodesOldToNew[j]; // new vtkId --> old SMDS id
351     }
352 }
353
354 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
355                                      std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
356                                      vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
357                                      int start, int end)
358 {
359   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
360   for (int j = start; j < end; j++)
361     {
362       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
363       idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
364       vtkIdType oldLoc = this->Locations->GetValue(j);
365       vtkIdType nbpts;
366       vtkIdType *oldPtsCell = 0;
367       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
368       assert(nbpts < NBMAXNODESINCELL);
369       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
370       for (int l = 0; l < nbpts; l++)
371         {
372           int oldval = oldPtsCell[l];
373           pointsCell[l] = idNodesOldToNew[oldval];
374           //MESSAGE("   " << oldval << " " << pointsCell[l]);
375         }
376       int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
377       int newLoc = newConnectivity->GetInsertLocation(nbpts);
378       //MESSAGE(newcnt << " " << newLoc);
379       newLocations->SetValue(alreadyCopied, newLoc);
380       alreadyCopied++;
381     }
382 }
383
384 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
385 {
386   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
387   return _cellIdToDownId[vtkCellId];
388 }
389
390 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
391 {
392   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
393   _cellIdToDownId[vtkCellId] = downId;
394 }
395
396 /*! Build downward connectivity: to do only when needed because heavy memory load.
397  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
398  *
399  */
400 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
401 {
402   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
403   // TODO calcul partiel sans edges
404
405   // --- erase previous data if any
406
407   for (int i = 0; i < _downArray.size(); i++)
408     {
409       if (_downArray[i])
410         delete _downArray[i];
411       _downArray[i] = 0;
412     }
413   _cellIdToDownId.clear();
414
415   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
416
417   _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
418
419   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
420   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
421   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
422   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
423   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
424   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
425   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
426   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
427   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
428   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
429   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
430   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
431   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
432   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
433
434   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
435
436   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
437
438   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
439   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
440   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
441   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
442   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
443   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
444   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
445   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
446
447   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
448   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
449   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
450   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
451   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
452   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
453
454   int GuessSize[VTK_QUADRATIC_TETRA];
455   GuessSize[VTK_LINE] = nbLineGuess;
456   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
457   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
458   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
459   GuessSize[VTK_QUAD] = nbLinQuadGuess;
460   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
461   GuessSize[VTK_TETRA] = nbLinTetra;
462   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
463   GuessSize[VTK_PYRAMID] = nbLinPyra;
464   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
465   GuessSize[VTK_WEDGE] = nbLinPrism;
466   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
467   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
468   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
469
470   _downArray[VTK_LINE]->allocate(nbLineGuess);
471   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
472   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
473   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
474   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
475   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
476   _downArray[VTK_TETRA]->allocate(nbLinTetra);
477   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
478   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
479   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
480   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
481   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
482   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
483   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
484
485   // --- iteration on vtkUnstructuredGrid cells, only faces
486   //     for each vtk face:
487   //       create a downward face entry with its downward id.
488   //       compute vtk volumes, create downward volumes entry.
489   //       mark face in downward volumes
490   //       mark volumes in downward face
491
492   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
493   int cellSize = this->Types->GetNumberOfTuples();
494   _cellIdToDownId.resize(cellSize, -1);
495
496   for (int i = 0; i < cellSize; i++)
497     {
498       int vtkFaceType = this->GetCellType(i);
499       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
500         {
501           int vtkFaceId = i;
502           //ASSERT(_downArray[vtkFaceType]);
503           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
504           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
505           downFace->setTempNodes(connFaceId, vtkFaceId);
506           int vols[2] = { -1, -1 };
507           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
508           //MESSAGE("nbVolumes="<< nbVolumes);
509           for (int ivol = 0; ivol < nbVolumes; ivol++)
510             {
511               int vtkVolId = vols[ivol];
512               int vtkVolType = this->GetCellType(vtkVolId);
513               //ASSERT(_downArray[vtkVolType]);
514               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
515               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
516               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
517               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
518             }
519         }
520     }
521
522   // --- iteration on vtkUnstructuredGrid cells, only volumes
523   //     for each vtk volume:
524   //       create downward volumes entry if not already done
525   //       build a temporary list of faces described with their nodes
526   //       for each face
527   //         compute the vtk volumes containing this face
528   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
529   //         if not, create a downward face entry (resizing of structure required)
530   //         (the downward faces store a temporary list of nodes to ease the comparison)
531   //         create downward volumes entry if not already done
532   //         mark volumes in downward face
533   //         mark face in downward volumes
534
535   CHRONOSTOP(20);
536   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
537
538   for (int i = 0; i < cellSize; i++)
539     {
540       int vtkType = this->GetCellType(i);
541       if (SMDS_Downward::getCellDimension(vtkType) == 3)
542         {
543           //CHRONO(31);
544           int vtkVolId = i;
545           // MESSAGE("vtk volume " << vtkVolId);
546           //ASSERT(_downArray[vtkType]);
547           int connVolId = _downArray[vtkType]->addCell(vtkVolId);
548
549           // --- find all the faces of the volume, describe the faces by their nodes
550
551           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
552           ListElemByNodesType facesWithNodes;
553           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
554           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
555           //CHRONOSTOP(31);
556           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
557             {
558               // --- find the volumes containing the face
559
560               //CHRONO(32);
561               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
562               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
563               int vols[2] = { -1, -1 };
564               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
565               int lg = facesWithNodes.elems[iface].nbNodes;
566               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
567               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
568
569               // --- check if face is registered in the volumes
570               //CHRONOSTOP(32);
571
572               //CHRONO(33);
573               int connFaceId = -1;
574               for (int ivol = 0; ivol < nbVolumes; ivol++)
575                 {
576                   int vtkVolId2 = vols[ivol];
577                   int vtkVolType = this->GetCellType(vtkVolId2);
578                   //ASSERT(_downArray[vtkVolType]);
579                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
580                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
581                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
582                   if (connFaceId >= 0)
583                     break; // --- face already created
584                 }//CHRONOSTOP(33);
585
586               // --- if face is not registered in the volumes, create face
587
588               //CHRONO(34);
589               if (connFaceId < 0)
590                 {
591                   connFaceId = _downArray[vtkFaceType]->addCell();
592                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
593                 }//CHRONOSTOP(34);
594
595               // --- mark volumes in downward face and mark face in downward volumes
596
597               //CHRONO(35);
598               for (int ivol = 0; ivol < nbVolumes; ivol++)
599                 {
600                   int vtkVolId2 = vols[ivol];
601                   int vtkVolType = this->GetCellType(vtkVolId2);
602                   //ASSERT(_downArray[vtkVolType]);
603                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
604                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
605                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
606                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
607                 }//CHRONOSTOP(35);
608             }
609         }
610     }
611
612   // --- iteration on vtkUnstructuredGrid cells, only edges
613   //     for each vtk edge:
614   //       create downward edge entry
615   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
616   //       find downward faces
617   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
618   //       mark edge in downward faces
619   //       mark faces in downward edge
620
621   CHRONOSTOP(21);
622   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
623
624   for (int i = 0; i < cellSize; i++)
625     {
626       int vtkEdgeType = this->GetCellType(i);
627       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
628         {
629           int vtkEdgeId = i;
630           //ASSERT(_downArray[vtkEdgeType]);
631           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
632           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
633           downEdge->setNodes(connEdgeId, vtkEdgeId);
634           vector<int> vtkIds;
635           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
636           int downFaces[1000];
637           unsigned char downTypes[1000];
638           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
639           for (int n = 0; n < nbDownFaces; n++)
640             {
641               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
642               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
643             }
644         }
645     }
646
647   // --- iteration on downward faces (they are all listed now)
648   //     for each downward face:
649   //       build a temporary list of edges with their ordered list of nodes
650   //       for each edge:
651   //         find all the vtk cells containing this edge
652   //         then identify all the downward faces containing the edge, from the vtk cells
653   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
654   //         if not, create a downward edge entry with the node id's
655   //         mark edge in downward faces
656   //         mark downward faces in edge (size of list unknown, to be allocated)
657
658   CHRONOSTOP(22);CHRONO(23);
659
660   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
661     {
662       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
663         continue;
664
665       // --- find all the edges of the face, describe the edges by their nodes
666
667       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
668       int maxId = downFace->getMaxId();
669       for (int faceId = 0; faceId < maxId; faceId++)
670         {
671           //CHRONO(40);
672           ListElemByNodesType edgesWithNodes;
673           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
674           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
675
676           //CHRONOSTOP(40);
677           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
678             {
679
680               // --- check if the edge is already registered by exploration of the faces
681
682               //CHRONO(41);
683               vector<int> vtkIds;
684               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
685               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
686               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
687               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
688               //CHRONOSTOP(41);CHRONO(42);
689               int downFaces[1000];
690               unsigned char downTypes[1000];
691               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
692               //CHRONOSTOP(42);
693
694               //CHRONO(43);
695               int connEdgeId = -1;
696               for (int idf = 0; idf < nbDownFaces; idf++)
697                 {
698                   int faceId2 = downFaces[idf];
699                   int faceType = downTypes[idf];
700                   //ASSERT(_downArray[faceType]);
701                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
702                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
703                   if (connEdgeId >= 0)
704                     break; // --- edge already created
705                 }//CHRONOSTOP(43);
706
707               // --- if edge is not registered in the faces, create edge
708
709               if (connEdgeId < 0)
710                 {
711                   //CHRONO(44);
712                   connEdgeId = _downArray[vtkEdgeType]->addCell();
713                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
714                   //CHRONOSTOP(44);
715                 }
716
717               // --- mark faces in downward edge and mark edge in downward faces
718
719               //CHRONO(45);
720               for (int idf = 0; idf < nbDownFaces; idf++)
721                 {
722                   int faceId2 = downFaces[idf];
723                   int faceType = downTypes[idf];
724                   //ASSERT(_downArray[faceType]);
725                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
726                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
727                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
728                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
729                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
730                 }//CHRONOSTOP(45);
731             }
732         }
733     }
734
735   CHRONOSTOP(23);CHRONO(24);
736
737   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
738   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
739
740   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
741     {
742       if (SMDS_Downward *down = _downArray[vtkType])
743         {
744           down->compactStorage();
745         }
746     }
747
748   // --- Statistics
749
750   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
751     {
752       if (SMDS_Downward *down = _downArray[vtkType])
753         {
754           if (down->getMaxId())
755             {
756               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
757                   << GuessSize[vtkType] << " real: " << down->getMaxId());
758             }
759         }
760     }CHRONOSTOP(24);CHRONOSTOP(2);
761   counters::stats();
762 }
763
764 /*! Get the neighbors of a cell.
765  * Only the neighbors having the dimension of the cell are taken into account
766  * (neighbors of a volume are the volumes sharing a face with this volume,
767  *  neighbors of a face are the faces sharing an edge with this face...).
768  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
769  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
770  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
771  * @param vtkId the vtk id of the cell
772  * @return number of neighbors
773  */
774 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
775 {
776   int vtkType = this->GetCellType(vtkId);
777   int cellDim = SMDS_Downward::getCellDimension(vtkType);
778   if (cellDim != 3)
779     return 0; // TODO voisins des faces ou edges
780   int cellId = this->_cellIdToDownId[vtkId];
781
782   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
783   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
784   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
785
786   // --- iteration on faces of the 3D cell.
787
788   int nb = 0;
789   for (int i = 0; i < nbCells; i++)
790     {
791       int downId = downCells[i];
792       int cellType = downTyp[i];
793       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
794       const int *upCells = _downArray[cellType]->getUpCells(downId);
795       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
796
797       // --- max 2 upCells, one is this cell, the other is a neighbor
798
799       for (int j = 0; j < nbUp; j++)
800         {
801           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
802             continue;
803           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
804           neighborsVtkIds[nb] = vtkNeighbor;
805           downIds[nb] = downId;
806           downTypes[nb] = cellType;
807           nb++;
808         }
809       if (nb >= NBMAXNEIGHBORS)
810         assert(0);
811     }
812   return nb;
813 }
814
815 /*! get the node id's of a cell.
816  * The cell is defined by it's downward connectivity id and type.
817  * @param nodeSet set of of vtk node id's to fill.
818  * @param downId downward connectivity id of the cell.
819  * @param downType type of cell.
820  */
821 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
822 {
823   _downArray[downType]->getNodeIds(downId, nodeSet);
824 }
825
826 /*! change some nodes in cell without modifying type or internal connectivity.
827  * Nodes inverse connectivity is maintained up to date.
828  * @param vtkVolId vtk id of the cell
829  * @param localClonedNodeIds map old node id to new node id.
830  */
831 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
832 {
833   vtkIdType npts = 0;
834   vtkIdType *pts; // will refer to the point id's of the face
835   this->GetCellPoints(vtkVolId, npts, pts);
836   for (int i = 0; i < npts; i++)
837     {
838       if (localClonedNodeIds.count(pts[i]))
839         {
840           vtkIdType oldpt = pts[i];
841           pts[i] = localClonedNodeIds[oldpt];
842           //MESSAGE(oldpt << " --> " << pts[i]);
843           //this->RemoveReferenceToCell(oldpt, vtkVolId);
844           //this->AddReferenceToCell(pts[i], vtkVolId);
845         }
846     }
847 }
848
849 /*! Create a volume (prism or hexahedron) by duplication of a face.
850  * the nodes of the new face are already created.
851  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
852  * @param localClonedNodeIds map old node id to new node id.
853  * @return vtk id of the new volume.
854  */
855 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<int>& orderedNodes)
856 {
857   int vtkType = this->GetCellType(vtkVolId);
858   int cellDim = SMDS_Downward::getCellDimension(vtkType);
859   if (cellDim != 3)
860     return 0;
861   SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
862   int downVolId = this->_cellIdToDownId[vtkVolId];
863   downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
864   return orderedNodes.size();
865 }
866
867 void SMDS_UnstructuredGrid::BuildLinks()
868 {
869   // Remove the old links if they are already built
870   if (this->Links)
871     {
872     this->Links->UnRegister(this);
873     }
874
875   this->Links = SMDS_CellLinks::New();
876   this->Links->Allocate(this->GetNumberOfPoints());
877   this->Links->Register(this);
878   this->Links->BuildLinks(this, this->Connectivity);
879   this->Links->Delete();
880 }