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