Salome HOME
Replace oe by ?
[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 #include "SMDS_MeshVolume.hxx"
7
8 #include "utilities.h"
9
10 #include <vtkCellArray.h>
11 #include <vtkCellLinks.h>
12 #include <vtkIdTypeArray.h>
13 #include <vtkUnsignedCharArray.h>
14
15 #include <list>
16 #include <climits>
17
18 using namespace std;
19
20 SMDS_CellLinks* SMDS_CellLinks::New()
21 {
22   MESSAGE("SMDS_CellLinks::New");
23   return new SMDS_CellLinks();
24 }
25
26 vtkCellLinks::Link* SMDS_CellLinks::ResizeL(vtkIdType sz)
27 {
28   return vtkCellLinks::Resize(sz);
29 }
30
31 vtkIdType SMDS_CellLinks::GetLinksSize()
32 {
33   return this->Size;
34 }
35
36 SMDS_CellLinks::SMDS_CellLinks() :
37   vtkCellLinks()
38 {
39 }
40
41 SMDS_CellLinks::~SMDS_CellLinks()
42 {
43 }
44
45 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
46 {
47   MESSAGE("SMDS_UnstructuredGrid::New");
48   return new SMDS_UnstructuredGrid();
49 }
50
51 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
52   vtkUnstructuredGrid()
53 {
54   _cellIdToDownId.clear();
55   _downTypes.clear();
56   _downArray.clear();
57   _mesh = 0;
58 }
59
60 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
61 {
62 }
63
64 unsigned long SMDS_UnstructuredGrid::GetMTime()
65 {
66   unsigned long mtime = vtkUnstructuredGrid::GetMTime();
67   MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
68   return mtime;
69 }
70
71 void SMDS_UnstructuredGrid::Update()
72 {
73   MESSAGE("SMDS_UnstructuredGrid::Update");
74   return vtkUnstructuredGrid::Update();
75 }
76
77 void SMDS_UnstructuredGrid::UpdateInformation()
78 {
79   MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
80   return vtkUnstructuredGrid::UpdateInformation();
81 }
82
83 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
84 {
85   // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
86   //MESSAGE("*********************** SMDS_UnstructuredGrid::GetPoints " << this->Points << " " << vtkUnstructuredGrid::GetPoints());
87   return this->Points;
88 }
89
90 //#ifdef VTK_HAVE_POLYHEDRON
91 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
92 {
93   if (type != VTK_POLYHEDRON)
94     return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
95
96   // --- type = VTK_POLYHEDRON
97   //MESSAGE("InsertNextLinkedCell VTK_POLYHEDRON");
98   int cellid = this->InsertNextCell(type, npts, pts);
99
100   set<vtkIdType> setOfNodes;
101   setOfNodes.clear();
102   int nbfaces = npts;
103   int i = 0;
104   for (int nf = 0; nf < nbfaces; nf++)
105     {
106       int nbnodes = pts[i];
107       i++;
108       for (int k = 0; k < nbnodes; k++)
109         {
110           //MESSAGE(" cell " << cellid << " face " << nf << " node " << pts[i]);
111           setOfNodes.insert(pts[i]);
112           i++;
113         }
114     }
115
116   set<vtkIdType>::iterator it = setOfNodes.begin();
117   for (; it != setOfNodes.end(); ++it)
118     {
119       //MESSAGE("reverse link for node " << *it << " cell " << cellid);
120       this->Links->ResizeCellList(*it, 1);
121       this->Links->AddCellReference(cellid, *it);
122     }
123
124   return cellid;
125 }
126 //#endif
127
128 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
129 {
130   _mesh = mesh;
131 }
132
133 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
134                                         std::vector<int>& idCellsOldToNew, int newCellSize)
135 {
136   MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
137   int alreadyCopied = 0;
138
139   // --- if newNodeSize, create a new compacted vtkPoints
140
141   vtkPoints *newPoints = vtkPoints::New();
142   newPoints->SetDataType(VTK_DOUBLE);
143   newPoints->SetNumberOfPoints(newNodeSize);
144   if (newNodeSize)
145     {
146       MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
147       // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
148       // using double type for storing coordinates of nodes instead float.
149       int oldNodeSize = idNodesOldToNew.size();
150
151       int i = 0;
152       while ( i < oldNodeSize )
153       {
154         // skip a hole if any
155         while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
156           ++i;
157         int startBloc = i;
158         // look for a block end
159         while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
160           ++i;
161         int endBloc = i;
162         copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
163       }
164       newPoints->Squeeze();
165     }
166
167   // --- create new compacted Connectivity, Locations and Types
168
169   int oldCellSize = this->Types->GetNumberOfTuples();
170
171   vtkCellArray *newConnectivity = vtkCellArray::New();
172   newConnectivity->Initialize();
173   int oldCellDataSize = this->Connectivity->GetData()->GetSize();
174   newConnectivity->Allocate(oldCellDataSize);
175   MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
176
177   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
178   newTypes->Initialize();
179   newTypes->SetNumberOfValues(newCellSize);
180
181   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
182   newLocations->Initialize();
183   newLocations->SetNumberOfValues(newCellSize);
184
185   // TODO some polyhedron may be huge (only in some tests)
186   vtkIdType tmpid[NBMAXNODESINCELL];
187   vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
188
189   alreadyCopied = 0;
190   int i = 0;
191   while ( i < oldCellSize )
192   {
193     // skip a hole if any
194     while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
195       ++i;
196     int startBloc = i;
197     // look for a block end
198     while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
199       ++i;
200     int endBloc = i;
201     if ( endBloc > startBloc )
202       copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
203                alreadyCopied, startBloc, endBloc);
204   }
205
206   newConnectivity->Squeeze();
207
208   if (1/*newNodeSize*/)
209     {
210       MESSAGE("------- newNodeSize, setPoints");
211       this->SetPoints(newPoints);
212       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
213     }
214
215   if (this->FaceLocations)
216     {
217       vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
218       newFaceLocations->Initialize();
219       newFaceLocations->Allocate(newTypes->GetSize());
220       vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
221       newFaces->Initialize();
222       newFaces->Allocate(this->Faces->GetSize());
223       for (int i = 0; i < oldCellSize; i++)
224         {
225           if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
226             continue;
227           int newCellId = idCellsOldToNew[i];
228           if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
229             {
230                newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
231                int oldFaceLoc = this->FaceLocations->GetValue(i);
232                int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
233                newFaces->InsertNextValue(nCellFaces);
234                for (int n=0; n<nCellFaces; n++)
235                  {
236                    int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
237                    newFaces->InsertNextValue(nptsInFace);
238                    for (int k=0; k<nptsInFace; k++)
239                      {
240                        int oldpt = this->Faces->GetValue(oldFaceLoc++);
241                        newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
242                      }
243                  }
244             }
245           else
246             {
247                newFaceLocations->InsertNextValue(-1);
248             }
249         }
250       newFaceLocations->Squeeze();
251       newFaces->Squeeze();
252       newFaceLocations->Register(this);
253       newFaces->Register(this);
254       this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
255       newFaceLocations->Delete();
256       newFaces->Delete();
257     }
258   else
259     this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
260
261   newPoints->Delete();
262   newTypes->Delete();
263   newLocations->Delete();
264   newConnectivity->Delete();
265   this->BuildLinks();
266 }
267
268 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
269                                       int start, int end)
270 {
271   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
272   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
273   void *source = this->Points->GetVoidPointer(3 * start);
274   int nbPoints = end - start;
275   if (nbPoints > 0)
276     {
277       memcpy(target, source, 3 * sizeof(double) * nbPoints);
278       for (int j = start; j < end; j++)
279         idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
280     }
281 }
282
283 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
284                                      std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
285                                      vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
286                                      int start, int end)
287 {
288   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
289   for (int j = start; j < end; j++)
290     {
291       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
292       idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
293       vtkIdType oldLoc = this->Locations->GetValue(j);
294       vtkIdType nbpts;
295       vtkIdType *oldPtsCell = 0;
296       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
297       assert(nbpts < NBMAXNODESINCELL);
298       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
299       for (int l = 0; l < nbpts; l++)
300         {
301           int oldval = oldPtsCell[l];
302           pointsCell[l] = idNodesOldToNew[oldval];
303           //MESSAGE("   " << oldval << " " << pointsCell[l]);
304         }
305       /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
306       int newLoc = newConnectivity->GetInsertLocation(nbpts);
307       //MESSAGE(newcnt << " " << newLoc);
308       newLocations->SetValue(alreadyCopied, newLoc);
309       alreadyCopied++;
310     }
311 }
312
313 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
314 {
315   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
316   return _cellIdToDownId[vtkCellId];
317 }
318
319 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
320 {
321   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
322   _cellIdToDownId[vtkCellId] = downId;
323 }
324
325 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
326 {
327   for (int i = 0; i < _downArray.size(); i++)
328     {
329       if (_downArray[i])
330         delete _downArray[i];
331       _downArray[i] = 0;
332     }
333   _cellIdToDownId.clear();
334 }
335
336 /*! Build downward connectivity: to do only when needed because heavy memory load.
337  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
338  *
339  */
340 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
341 {
342   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
343   // TODO calcul partiel sans edges
344
345   // --- erase previous data if any
346
347   this->CleanDownwardConnectivity();
348
349   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
350
351   _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
352
353   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
354   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
355   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
356   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
357   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
358   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
359   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
360   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
361   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
362   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
363   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
364   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
365   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
366   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
367
368   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
369
370   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
371
372   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
373   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
374   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
375   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
376   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
377   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
378   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
379   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
380
381   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
382   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
383   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
384   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
385   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
386   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
387
388   int GuessSize[VTK_QUADRATIC_TETRA];
389   GuessSize[VTK_LINE] = nbLineGuess;
390   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
391   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
392   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
393   GuessSize[VTK_QUAD] = nbLinQuadGuess;
394   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
395   GuessSize[VTK_TETRA] = nbLinTetra;
396   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
397   GuessSize[VTK_PYRAMID] = nbLinPyra;
398   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
399   GuessSize[VTK_WEDGE] = nbLinPrism;
400   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
401   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
402   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
403
404   _downArray[VTK_LINE]->allocate(nbLineGuess);
405   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
406   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
407   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
408   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
409   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
410   _downArray[VTK_TETRA]->allocate(nbLinTetra);
411   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
412   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
413   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
414   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
415   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
416   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
417   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
418
419   // --- iteration on vtkUnstructuredGrid cells, only faces
420   //     for each vtk face:
421   //       create a downward face entry with its downward id.
422   //       compute vtk volumes, create downward volumes entry.
423   //       mark face in downward volumes
424   //       mark volumes in downward face
425
426   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
427   int cellSize = this->Types->GetNumberOfTuples();
428   _cellIdToDownId.resize(cellSize, -1);
429
430   for (int i = 0; i < cellSize; i++)
431     {
432       int vtkFaceType = this->GetCellType(i);
433       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
434         {
435           int vtkFaceId = i;
436           //ASSERT(_downArray[vtkFaceType]);
437           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
438           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
439           downFace->setTempNodes(connFaceId, vtkFaceId);
440           int vols[2] = { -1, -1 };
441           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
442           //MESSAGE("nbVolumes="<< nbVolumes);
443           for (int ivol = 0; ivol < nbVolumes; ivol++)
444             {
445               int vtkVolId = vols[ivol];
446               int vtkVolType = this->GetCellType(vtkVolId);
447               //ASSERT(_downArray[vtkVolType]);
448               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
449               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
450               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
451               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
452             }
453         }
454     }
455
456   // --- iteration on vtkUnstructuredGrid cells, only volumes
457   //     for each vtk volume:
458   //       create downward volumes entry if not already done
459   //       build a temporary list of faces described with their nodes
460   //       for each face
461   //         compute the vtk volumes containing this face
462   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
463   //         if not, create a downward face entry (resizing of structure required)
464   //         (the downward faces store a temporary list of nodes to ease the comparison)
465   //         create downward volumes entry if not already done
466   //         mark volumes in downward face
467   //         mark face in downward volumes
468
469   CHRONOSTOP(20);
470   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
471
472   for (int i = 0; i < cellSize; i++)
473     {
474       int vtkType = this->GetCellType(i);
475       if (SMDS_Downward::getCellDimension(vtkType) == 3)
476         {
477           //CHRONO(31);
478           int vtkVolId = i;
479           // MESSAGE("vtk volume " << vtkVolId);
480           //ASSERT(_downArray[vtkType]);
481           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
482
483           // --- find all the faces of the volume, describe the faces by their nodes
484
485           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
486           ListElemByNodesType facesWithNodes;
487           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
488           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
489           //CHRONOSTOP(31);
490           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
491             {
492               // --- find the volumes containing the face
493
494               //CHRONO(32);
495               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
496               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
497               int vols[2] = { -1, -1 };
498               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
499               int lg = facesWithNodes.elems[iface].nbNodes;
500               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
501               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
502
503               // --- check if face is registered in the volumes
504               //CHRONOSTOP(32);
505
506               //CHRONO(33);
507               int connFaceId = -1;
508               for (int ivol = 0; ivol < nbVolumes; ivol++)
509                 {
510                   int vtkVolId2 = vols[ivol];
511                   int vtkVolType = this->GetCellType(vtkVolId2);
512                   //ASSERT(_downArray[vtkVolType]);
513                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
514                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
515                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
516                   if (connFaceId >= 0)
517                     break; // --- face already created
518                 }//CHRONOSTOP(33);
519
520               // --- if face is not registered in the volumes, create face
521
522               //CHRONO(34);
523               if (connFaceId < 0)
524                 {
525                   connFaceId = _downArray[vtkFaceType]->addCell();
526                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
527                 }//CHRONOSTOP(34);
528
529               // --- mark volumes in downward face and mark face in downward volumes
530
531               //CHRONO(35);
532               for (int ivol = 0; ivol < nbVolumes; ivol++)
533                 {
534                   int vtkVolId2 = vols[ivol];
535                   int vtkVolType = this->GetCellType(vtkVolId2);
536                   //ASSERT(_downArray[vtkVolType]);
537                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
538                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
539                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
540                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
541                 }//CHRONOSTOP(35);
542             }
543         }
544     }
545
546   // --- iteration on vtkUnstructuredGrid cells, only edges
547   //     for each vtk edge:
548   //       create downward edge entry
549   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
550   //       find downward faces
551   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
552   //       mark edge in downward faces
553   //       mark faces in downward edge
554
555   CHRONOSTOP(21);
556   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
557
558   for (int i = 0; i < cellSize; i++)
559     {
560       int vtkEdgeType = this->GetCellType(i);
561       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
562         {
563           int vtkEdgeId = i;
564           //ASSERT(_downArray[vtkEdgeType]);
565           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
566           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
567           downEdge->setNodes(connEdgeId, vtkEdgeId);
568           vector<int> vtkIds;
569           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
570           int downFaces[1000];
571           unsigned char downTypes[1000];
572           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
573           for (int n = 0; n < nbDownFaces; n++)
574             {
575               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
576               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
577             }
578         }
579     }
580
581   // --- iteration on downward faces (they are all listed now)
582   //     for each downward face:
583   //       build a temporary list of edges with their ordered list of nodes
584   //       for each edge:
585   //         find all the vtk cells containing this edge
586   //         then identify all the downward faces containing the edge, from the vtk cells
587   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
588   //         if not, create a downward edge entry with the node id's
589   //         mark edge in downward faces
590   //         mark downward faces in edge (size of list unknown, to be allocated)
591
592   CHRONOSTOP(22);CHRONO(23);
593
594   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
595     {
596       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
597         continue;
598
599       // --- find all the edges of the face, describe the edges by their nodes
600
601       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
602       int maxId = downFace->getMaxId();
603       for (int faceId = 0; faceId < maxId; faceId++)
604         {
605           //CHRONO(40);
606           ListElemByNodesType edgesWithNodes;
607           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
608           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
609
610           //CHRONOSTOP(40);
611           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
612             {
613
614               // --- check if the edge is already registered by exploration of the faces
615
616               //CHRONO(41);
617               vector<int> vtkIds;
618               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
619               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
620               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
621               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
622               //CHRONOSTOP(41);CHRONO(42);
623               int downFaces[1000];
624               unsigned char downTypes[1000];
625               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
626               //CHRONOSTOP(42);
627
628               //CHRONO(43);
629               int connEdgeId = -1;
630               for (int idf = 0; idf < nbDownFaces; idf++)
631                 {
632                   int faceId2 = downFaces[idf];
633                   int faceType = downTypes[idf];
634                   //ASSERT(_downArray[faceType]);
635                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
636                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
637                   if (connEdgeId >= 0)
638                     break; // --- edge already created
639                 }//CHRONOSTOP(43);
640
641               // --- if edge is not registered in the faces, create edge
642
643               if (connEdgeId < 0)
644                 {
645                   //CHRONO(44);
646                   connEdgeId = _downArray[vtkEdgeType]->addCell();
647                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
648                   //CHRONOSTOP(44);
649                 }
650
651               // --- mark faces in downward edge and mark edge in downward faces
652
653               //CHRONO(45);
654               for (int idf = 0; idf < nbDownFaces; idf++)
655                 {
656                   int faceId2 = downFaces[idf];
657                   int faceType = downTypes[idf];
658                   //ASSERT(_downArray[faceType]);
659                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
660                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
661                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
662                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
663                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
664                 }//CHRONOSTOP(45);
665             }
666         }
667     }
668
669   CHRONOSTOP(23);CHRONO(24);
670
671   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
672   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
673
674   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
675     {
676       if (SMDS_Downward *down = _downArray[vtkType])
677         {
678           down->compactStorage();
679         }
680     }
681
682   // --- Statistics
683
684   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
685     {
686       if (SMDS_Downward *down = _downArray[vtkType])
687         {
688           if (down->getMaxId())
689             {
690               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
691                   << GuessSize[vtkType] << " real: " << down->getMaxId());
692             }
693         }
694     }CHRONOSTOP(24);CHRONOSTOP(2);
695   counters::stats();
696 }
697
698 /*! Get the neighbors of a cell.
699  * Only the neighbors having the dimension of the cell are taken into account
700  * (neighbors of a volume are the volumes sharing a face with this volume,
701  *  neighbors of a face are the faces sharing an edge with this face...).
702  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
703  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
704  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
705  * @param vtkId the vtk id of the cell
706  * @return number of neighbors
707  */
708 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
709 {
710   int vtkType = this->GetCellType(vtkId);
711   int cellDim = SMDS_Downward::getCellDimension(vtkType);
712   if (cellDim <2)
713     return 0; // TODO voisins des edges = edges connectees
714   int cellId = this->_cellIdToDownId[vtkId];
715
716   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
717   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
718   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
719
720   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
721
722   int nb = 0;
723   for (int i = 0; i < nbCells; i++)
724     {
725       int downId = downCells[i];
726       int cellType = downTyp[i];
727       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
728       const int *upCells = _downArray[cellType]->getUpCells(downId);
729       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
730
731       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
732       //    for a face, number of neighbors (connected faces) not known
733
734       for (int j = 0; j < nbUp; j++)
735         {
736           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
737             continue;
738           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
739           neighborsVtkIds[nb] = vtkNeighbor;
740           downIds[nb] = downId;
741           downTypes[nb] = cellType;
742           nb++;
743         }
744       if (nb >= NBMAXNEIGHBORS)
745         assert(0);
746     }
747   return nb;
748 }
749
750 /*! get the volumes containing a face or an edge of the grid
751  * The edge or face belongs to the vtkUnstructuredGrid
752  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
753  * @param vtkId vtk id of the face or edge
754  */
755 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
756 {
757   int vtkType = this->GetCellType(vtkId);
758   int dim = SMDS_Downward::getCellDimension(vtkType);
759   int nbFaces = 0;
760   int faces[1000];
761   unsigned char cellTypes[1000];
762   int downCellId[1000];
763   if (dim == 1)
764     {
765       int downId = this->CellIdToDownId(vtkId);
766       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
767       const int *upCells = _downArray[vtkType]->getUpCells(downId);
768       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
769       for (int i=0; i< nbFaces; i++)
770         {
771           faces[i] = _downArray[upTypes[i]]->getVtkCellId(upCells[i]);
772           cellTypes[i] = upTypes[i];
773           downCellId[i] = upCells[i];
774         }
775     }
776   else if (dim == 2)
777     {
778       nbFaces = 1;
779       faces[0] = vtkId;
780       cellTypes[0] = this->GetCellType(vtkId);
781       downCellId[0] = this->CellIdToDownId(vtkId);
782     }
783
784   int nbvol =0;
785   for (int i=0; i<nbFaces; i++)
786     {
787       int vtkTypeFace = cellTypes[i];
788       int downId = downCellId[i];
789       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
790       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
791       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
792        for (int j=0; j<nv; j++)
793         {
794           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
795           if (vtkVolId >= 0)
796             volVtkIds[nbvol++] = vtkVolId;
797         }
798     }
799   return nbvol;
800 }
801
802 /*! get the volumes containing a face or an edge of the downward structure
803  * The edge or face does not necessary belong to the vtkUnstructuredGrid
804  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
805  * @param downId id in the downward structure
806  * @param downType type of cell
807  */
808 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
809 {
810   int vtkType = downType;
811   int dim = SMDS_Downward::getCellDimension(vtkType);
812   int nbFaces = 0;
813   int faces[1000];
814   unsigned char cellTypes[1000];
815   int downCellId[1000];
816   if (dim == 1)
817     {
818       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
819       const int *upCells = _downArray[vtkType]->getUpCells(downId);
820       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
821       for (int i=0; i< nbFaces; i++)
822         {
823           faces[i] = _downArray[upTypes[i]]->getVtkCellId(upCells[i]);
824           cellTypes[i] = upTypes[i];
825           downCellId[i] = upCells[i];
826         }
827     }
828   else if (dim == 2)
829     {
830       nbFaces = 1;
831       cellTypes[0] = vtkType;
832       downCellId[0] = downId;
833     }
834
835   int nbvol =0;
836   for (int i=0; i<nbFaces; i++)
837     {
838       int vtkTypeFace = cellTypes[i];
839       int downId = downCellId[i];
840       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
841       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
842       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
843        for (int j=0; j<nv; j++)
844         {
845           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
846           if (vtkVolId >= 0)
847             volVtkIds[nbvol++] = vtkVolId;
848         }
849     }
850   return nbvol;
851 }
852
853 /*! get the node id's of a cell.
854  * The cell is defined by it's downward connectivity id and type.
855  * @param nodeSet set of of vtk node id's to fill.
856  * @param downId downward connectivity id of the cell.
857  * @param downType type of cell.
858  */
859 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
860 {
861   _downArray[downType]->getNodeIds(downId, nodeSet);
862 }
863
864 /*! change some nodes in cell without modifying type or internal connectivity.
865  * Nodes inverse connectivity is maintained up to date.
866  * @param vtkVolId vtk id of the cell
867  * @param localClonedNodeIds map old node id to new node id.
868  */
869 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
870 {
871   vtkIdType npts = 0;
872   vtkIdType *pts; // will refer to the point id's of the face
873   this->GetCellPoints(vtkVolId, npts, pts);
874   for (int i = 0; i < npts; i++)
875     {
876       if (localClonedNodeIds.count(pts[i]))
877         {
878           vtkIdType oldpt = pts[i];
879           pts[i] = localClonedNodeIds[oldpt];
880           //MESSAGE(oldpt << " --> " << pts[i]);
881           //this->RemoveReferenceToCell(oldpt, vtkVolId);
882           //this->AddReferenceToCell(pts[i], vtkVolId);
883         }
884     }
885 }
886
887 /*! reorder the nodes of a face
888  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
889  * @param orderedNodes list of nodes to reorder (in out)
890  * @return size of the list
891  */
892 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes)
893 {
894   int vtkType = this->GetCellType(vtkVolId);
895   int cellDim = SMDS_Downward::getCellDimension(vtkType);
896   if (cellDim != 3)
897     return 0;
898   SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
899   int downVolId = this->_cellIdToDownId[vtkVolId];
900   downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
901   return orderedNodes.size();
902 }
903
904 void SMDS_UnstructuredGrid::BuildLinks()
905 {
906   // Remove the old links if they are already built
907   if (this->Links)
908     {
909     this->Links->UnRegister(this);
910     }
911
912   this->Links = SMDS_CellLinks::New();
913   this->Links->Allocate(this->GetNumberOfPoints());
914   this->Links->Register(this);
915   this->Links->BuildLinks(this, this->Connectivity);
916   this->Links->Delete();
917 }
918
919 /*! Create a volume (prism or hexahedron) by duplication of a face.
920  * Designed for use in creation of flat elements separating volume domains.
921  * A face separating two domains is shared by two volume cells.
922  * All the nodes are already created (for the two faces).
923  * Each original Node is associated to corresponding nodes in the domains.
924  * Some nodes may be duplicated for more than two domains, when domain separations intersect.
925  * In that case, even some of the nodes to use for the original face may be changed.
926  * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
927  * @param domain1: domain of the original face
928  * @param domain2: domain of the duplicated face
929  * @param originalNodes: the vtk node ids of the original face
930  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
931  * @return ok if success.
932  */
933 SMDS_MeshVolume* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
934                                                   int domain1,
935                                                   int domain2,
936                                                   std::set<int>& originalNodes,
937                                                   std::map<int, std::map<int, int> >& nodeDomains,
938                                                   std::map<int, std::map<long, int> >& nodeQuadDomains)
939 {
940   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
941   vector<vtkIdType> orderedOriginals;
942   orderedOriginals.clear();
943   set<int>::const_iterator it = originalNodes.begin();
944   for (; it != originalNodes.end(); ++it)
945     orderedOriginals.push_back(*it);
946
947   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, orderedOriginals);
948   vector<vtkIdType> orderedNodes;
949
950   switch (orderedOriginals.size())
951   {
952     case 3:
953     case 4:
954       for (int i = 0; i < nbNodes; i++)
955         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
956       for (int i = 0; i < nbNodes; i++)
957         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
958       break;
959     case 6:
960     case 8:
961       {
962         long dom1 = domain1;
963         long dom2 = domain2;
964         long dom1_2; // for nodeQuadDomains
965         if (domain1 < domain2)
966           dom1_2 = dom1 + INT_MAX * dom2;
967         else
968           dom1_2 = dom2 + INT_MAX * dom1;
969         //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
970         int ima = orderedOriginals.size();
971         int mid = orderedOriginals.size() / 2;
972         //cerr << "ima=" << ima << " mid=" << mid << endl;
973         for (int i = 0; i < mid; i++)
974           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
975         for (int i = 0; i < mid; i++)
976           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
977         for (int i = mid; i < ima; i++)
978           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
979         for (int i = mid; i < ima; i++)
980           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
981         for (int i = 0; i < mid; i++)
982           {
983             int oldId = orderedOriginals[i];
984             int newId;
985             if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
986               newId = nodeQuadDomains[oldId][dom1_2];
987             else
988               {
989                 double *coords = this->GetPoint(oldId);
990                 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
991                 newId = newNode->getVtkId();
992                 std::map<long, int> emptyMap;
993                 nodeQuadDomains[oldId] = emptyMap;
994                 nodeQuadDomains[oldId][dom1_2] = newId;
995               }
996             orderedNodes.push_back(newId);
997           }
998       }
999       break;
1000     default:
1001       ASSERT(0);
1002   }
1003
1004   SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1005
1006   // TODO update subshape list of elements and nodes
1007   return vol;
1008 }