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