Salome HOME
961ca8be2e9ee05acf2b2c78f71957710a1e8f85
[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   if((vtkCellId < 0) || (vtkCellId >= _cellIdToDownId.size()))
316     {
317       //MESSAGE("SMDS_UnstructuredGrid::CellIdToDownId structure not up to date: vtkCellId="
318       //    << vtkCellId << " max="<< _cellIdToDownId.size());
319       return -1;
320     }
321   return _cellIdToDownId[vtkCellId];
322 }
323
324 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
325 {
326   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
327   _cellIdToDownId[vtkCellId] = downId;
328 }
329
330 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
331 {
332   for (int i = 0; i < _downArray.size(); i++)
333     {
334       if (_downArray[i])
335         delete _downArray[i];
336       _downArray[i] = 0;
337     }
338   _cellIdToDownId.clear();
339 }
340
341 /*! Build downward connectivity: to do only when needed because heavy memory load.
342  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
343  *
344  */
345 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
346 {
347   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
348   // TODO calcul partiel sans edges
349
350   // --- erase previous data if any
351
352   this->CleanDownwardConnectivity();
353
354   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
355
356   _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
357
358   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
359   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
360   _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
361   _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
362   _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
363   _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
364   _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
365   _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
366   _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
367   _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
368   _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
369   _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
370   _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
371   _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
372
373   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
374
375   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
376
377   int nbLinTetra = meshInfo.NbTetras(ORDER_LINEAR);
378   int nbQuadTetra = meshInfo.NbTetras(ORDER_QUADRATIC);
379   int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
380   int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
381   int nbLinPrism = meshInfo.NbPrisms(ORDER_LINEAR);
382   int nbQuadPrism = meshInfo.NbPrisms(ORDER_QUADRATIC);
383   int nbLinHexa = meshInfo.NbHexas(ORDER_LINEAR);
384   int nbQuadHexa = meshInfo.NbHexas(ORDER_QUADRATIC);
385
386   int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
387   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
388   int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
389   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
390   int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
391   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
392
393   int GuessSize[VTK_QUADRATIC_TETRA];
394   GuessSize[VTK_LINE] = nbLineGuess;
395   GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
396   GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
397   GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
398   GuessSize[VTK_QUAD] = nbLinQuadGuess;
399   GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
400   GuessSize[VTK_TETRA] = nbLinTetra;
401   GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
402   GuessSize[VTK_PYRAMID] = nbLinPyra;
403   GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
404   GuessSize[VTK_WEDGE] = nbLinPrism;
405   GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
406   GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
407   GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
408
409   _downArray[VTK_LINE]->allocate(nbLineGuess);
410   _downArray[VTK_QUADRATIC_EDGE]->allocate(nbQuadEdgeGuess);
411   _downArray[VTK_TRIANGLE]->allocate(nbLinTriaGuess);
412   _downArray[VTK_QUADRATIC_TRIANGLE]->allocate(nbQuadTriaGuess);
413   _downArray[VTK_QUAD]->allocate(nbLinQuadGuess);
414   _downArray[VTK_QUADRATIC_QUAD]->allocate(nbQuadQuadGuess);
415   _downArray[VTK_TETRA]->allocate(nbLinTetra);
416   _downArray[VTK_QUADRATIC_TETRA]->allocate(nbQuadTetra);
417   _downArray[VTK_PYRAMID]->allocate(nbLinPyra);
418   _downArray[VTK_QUADRATIC_PYRAMID]->allocate(nbQuadPyra);
419   _downArray[VTK_WEDGE]->allocate(nbLinPrism);
420   _downArray[VTK_QUADRATIC_WEDGE]->allocate(nbQuadPrism);
421   _downArray[VTK_HEXAHEDRON]->allocate(nbLinHexa);
422   _downArray[VTK_QUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
423
424   // --- iteration on vtkUnstructuredGrid cells, only faces
425   //     for each vtk face:
426   //       create a downward face entry with its downward id.
427   //       compute vtk volumes, create downward volumes entry.
428   //       mark face in downward volumes
429   //       mark volumes in downward face
430
431   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
432   int cellSize = this->Types->GetNumberOfTuples();
433   _cellIdToDownId.resize(cellSize, -1);
434
435   for (int i = 0; i < cellSize; i++)
436     {
437       int vtkFaceType = this->GetCellType(i);
438       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
439         {
440           int vtkFaceId = i;
441           //ASSERT(_downArray[vtkFaceType]);
442           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
443           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
444           downFace->setTempNodes(connFaceId, vtkFaceId);
445           int vols[2] = { -1, -1 };
446           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
447           //MESSAGE("nbVolumes="<< nbVolumes);
448           for (int ivol = 0; ivol < nbVolumes; ivol++)
449             {
450               int vtkVolId = vols[ivol];
451               int vtkVolType = this->GetCellType(vtkVolId);
452               //ASSERT(_downArray[vtkVolType]);
453               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
454               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
455               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
456               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
457             }
458         }
459     }
460
461   // --- iteration on vtkUnstructuredGrid cells, only volumes
462   //     for each vtk volume:
463   //       create downward volumes entry if not already done
464   //       build a temporary list of faces described with their nodes
465   //       for each face
466   //         compute the vtk volumes containing this face
467   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
468   //         if not, create a downward face entry (resizing of structure required)
469   //         (the downward faces store a temporary list of nodes to ease the comparison)
470   //         create downward volumes entry if not already done
471   //         mark volumes in downward face
472   //         mark face in downward volumes
473
474   CHRONOSTOP(20);
475   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
476
477   for (int i = 0; i < cellSize; i++)
478     {
479       int vtkType = this->GetCellType(i);
480       if (SMDS_Downward::getCellDimension(vtkType) == 3)
481         {
482           //CHRONO(31);
483           int vtkVolId = i;
484           // MESSAGE("vtk volume " << vtkVolId);
485           //ASSERT(_downArray[vtkType]);
486           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
487
488           // --- find all the faces of the volume, describe the faces by their nodes
489
490           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
491           ListElemByNodesType facesWithNodes;
492           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
493           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
494           //CHRONOSTOP(31);
495           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
496             {
497               // --- find the volumes containing the face
498
499               //CHRONO(32);
500               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
501               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
502               int vols[2] = { -1, -1 };
503               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
504               int lg = facesWithNodes.elems[iface].nbNodes;
505               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
506               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
507
508               // --- check if face is registered in the volumes
509               //CHRONOSTOP(32);
510
511               //CHRONO(33);
512               int connFaceId = -1;
513               for (int ivol = 0; ivol < nbVolumes; ivol++)
514                 {
515                   int vtkVolId2 = vols[ivol];
516                   int vtkVolType = this->GetCellType(vtkVolId2);
517                   //ASSERT(_downArray[vtkVolType]);
518                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
519                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
520                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
521                   if (connFaceId >= 0)
522                     break; // --- face already created
523                 }//CHRONOSTOP(33);
524
525               // --- if face is not registered in the volumes, create face
526
527               //CHRONO(34);
528               if (connFaceId < 0)
529                 {
530                   connFaceId = _downArray[vtkFaceType]->addCell();
531                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
532                 }//CHRONOSTOP(34);
533
534               // --- mark volumes in downward face and mark face in downward volumes
535
536               //CHRONO(35);
537               for (int ivol = 0; ivol < nbVolumes; ivol++)
538                 {
539                   int vtkVolId2 = vols[ivol];
540                   int vtkVolType = this->GetCellType(vtkVolId2);
541                   //ASSERT(_downArray[vtkVolType]);
542                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
543                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
544                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
545                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
546                 }//CHRONOSTOP(35);
547             }
548         }
549     }
550
551   // --- iteration on vtkUnstructuredGrid cells, only edges
552   //     for each vtk edge:
553   //       create downward edge entry
554   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
555   //       find downward faces
556   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
557   //       mark edge in downward faces
558   //       mark faces in downward edge
559
560   CHRONOSTOP(21);
561   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
562
563   for (int i = 0; i < cellSize; i++)
564     {
565       int vtkEdgeType = this->GetCellType(i);
566       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
567         {
568           int vtkEdgeId = i;
569           //ASSERT(_downArray[vtkEdgeType]);
570           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
571           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
572           downEdge->setNodes(connEdgeId, vtkEdgeId);
573           vector<int> vtkIds;
574           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
575           int downFaces[1000];
576           unsigned char downTypes[1000];
577           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
578           for (int n = 0; n < nbDownFaces; n++)
579             {
580               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
581               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
582             }
583         }
584     }
585
586   // --- iteration on downward faces (they are all listed now)
587   //     for each downward face:
588   //       build a temporary list of edges with their ordered list of nodes
589   //       for each edge:
590   //         find all the vtk cells containing this edge
591   //         then identify all the downward faces containing the edge, from the vtk cells
592   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
593   //         if not, create a downward edge entry with the node id's
594   //         mark edge in downward faces
595   //         mark downward faces in edge (size of list unknown, to be allocated)
596
597   CHRONOSTOP(22);CHRONO(23);
598
599   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
600     {
601       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
602         continue;
603
604       // --- find all the edges of the face, describe the edges by their nodes
605
606       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
607       int maxId = downFace->getMaxId();
608       for (int faceId = 0; faceId < maxId; faceId++)
609         {
610           //CHRONO(40);
611           ListElemByNodesType edgesWithNodes;
612           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
613           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
614
615           //CHRONOSTOP(40);
616           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
617             {
618
619               // --- check if the edge is already registered by exploration of the faces
620
621               //CHRONO(41);
622               vector<int> vtkIds;
623               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
624               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
625               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
626               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
627               //CHRONOSTOP(41);CHRONO(42);
628               int downFaces[1000];
629               unsigned char downTypes[1000];
630               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
631               //CHRONOSTOP(42);
632
633               //CHRONO(43);
634               int connEdgeId = -1;
635               for (int idf = 0; idf < nbDownFaces; idf++)
636                 {
637                   int faceId2 = downFaces[idf];
638                   int faceType = downTypes[idf];
639                   //ASSERT(_downArray[faceType]);
640                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
641                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
642                   if (connEdgeId >= 0)
643                     break; // --- edge already created
644                 }//CHRONOSTOP(43);
645
646               // --- if edge is not registered in the faces, create edge
647
648               if (connEdgeId < 0)
649                 {
650                   //CHRONO(44);
651                   connEdgeId = _downArray[vtkEdgeType]->addCell();
652                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
653                   //CHRONOSTOP(44);
654                 }
655
656               // --- mark faces in downward edge and mark edge in downward faces
657
658               //CHRONO(45);
659               for (int idf = 0; idf < nbDownFaces; idf++)
660                 {
661                   int faceId2 = downFaces[idf];
662                   int faceType = downTypes[idf];
663                   //ASSERT(_downArray[faceType]);
664                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
665                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
666                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
667                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
668                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
669                 }//CHRONOSTOP(45);
670             }
671         }
672     }
673
674   CHRONOSTOP(23);CHRONO(24);
675
676   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
677   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
678
679   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
680     {
681       if (SMDS_Downward *down = _downArray[vtkType])
682         {
683           down->compactStorage();
684         }
685     }
686
687   // --- Statistics
688
689   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
690     {
691       if (SMDS_Downward *down = _downArray[vtkType])
692         {
693           if (down->getMaxId())
694             {
695               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
696                   << GuessSize[vtkType] << " real: " << down->getMaxId());
697             }
698         }
699     }CHRONOSTOP(24);CHRONOSTOP(2);
700   counters::stats();
701 }
702
703 /*! Get the neighbors of a cell.
704  * Only the neighbors having the dimension of the cell are taken into account
705  * (neighbors of a volume are the volumes sharing a face with this volume,
706  *  neighbors of a face are the faces sharing an edge with this face...).
707  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
708  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
709  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
710  * @param vtkId the vtk id of the cell
711  * @return number of neighbors
712  */
713 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
714 {
715   int vtkType = this->GetCellType(vtkId);
716   int cellDim = SMDS_Downward::getCellDimension(vtkType);
717   if (cellDim <2)
718     return 0; // TODO voisins des edges = edges connectees
719   int cellId = this->_cellIdToDownId[vtkId];
720
721   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
722   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
723   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
724
725   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
726
727   int nb = 0;
728   for (int i = 0; i < nbCells; i++)
729     {
730       int downId = downCells[i];
731       int cellType = downTyp[i];
732       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
733       const int *upCells = _downArray[cellType]->getUpCells(downId);
734       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
735
736       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
737       //    for a face, number of neighbors (connected faces) not known
738
739       for (int j = 0; j < nbUp; j++)
740         {
741           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
742             continue;
743           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
744           neighborsVtkIds[nb] = vtkNeighbor;
745           downIds[nb] = downId;
746           downTypes[nb] = cellType;
747           nb++;
748         }
749       if (nb >= NBMAXNEIGHBORS)
750         assert(0);
751     }
752   return nb;
753 }
754
755 /*! get the volumes containing a face or an edge of the grid
756  * The edge or face belongs to the vtkUnstructuredGrid
757  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
758  * @param vtkId vtk id of the face or edge
759  */
760 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
761 {
762   int vtkType = this->GetCellType(vtkId);
763   int dim = SMDS_Downward::getCellDimension(vtkType);
764   int nbFaces = 0;
765   unsigned char cellTypes[1000];
766   int downCellId[1000];
767   if (dim == 1)
768     {
769       int downId = this->CellIdToDownId(vtkId);
770       if (downId < 0)
771         {
772           MESSAGE("Downward structure not up to date: new edge not taken into account");
773           return 0;
774         }
775       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
776       const int *upCells = _downArray[vtkType]->getUpCells(downId);
777       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
778       for (int i=0; i< nbFaces; i++)
779         {
780           cellTypes[i] = upTypes[i];
781           downCellId[i] = upCells[i];
782         }
783     }
784   else if (dim == 2)
785     {
786       nbFaces = 1;
787       cellTypes[0] = this->GetCellType(vtkId);
788       int downId = this->CellIdToDownId(vtkId);
789       if (downId < 0)
790         {
791           MESSAGE("Downward structure not up to date: new face not taken into account");
792           return 0;
793         }
794       downCellId[0] = downId;
795     }
796
797   int nbvol =0;
798   for (int i=0; i<nbFaces; i++)
799     {
800       int vtkTypeFace = cellTypes[i];
801       int downId = downCellId[i];
802       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
803       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
804       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
805        for (int j=0; j<nv; j++)
806         {
807           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
808           if (vtkVolId >= 0)
809             volVtkIds[nbvol++] = vtkVolId;
810         }
811     }
812   return nbvol;
813 }
814
815 /*! get the volumes containing a face or an edge of the downward structure
816  * The edge or face does not necessary belong to the vtkUnstructuredGrid
817  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
818  * @param downId id in the downward structure
819  * @param downType type of cell
820  */
821 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
822 {
823   int vtkType = downType;
824   int dim = SMDS_Downward::getCellDimension(vtkType);
825   int nbFaces = 0;
826   unsigned char cellTypes[1000];
827   int downCellId[1000];
828   if (dim == 1)
829     {
830       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
831       const int *upCells = _downArray[vtkType]->getUpCells(downId);
832       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
833       for (int i=0; i< nbFaces; i++)
834         {
835           cellTypes[i] = upTypes[i];
836           downCellId[i] = upCells[i];
837         }
838     }
839   else if (dim == 2)
840     {
841       nbFaces = 1;
842       cellTypes[0] = vtkType;
843       downCellId[0] = downId;
844     }
845
846   int nbvol =0;
847   for (int i=0; i<nbFaces; i++)
848     {
849       int vtkTypeFace = cellTypes[i];
850       int downId = downCellId[i];
851       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
852       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
853       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
854        for (int j=0; j<nv; j++)
855         {
856           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
857           if (vtkVolId >= 0)
858             volVtkIds[nbvol++] = vtkVolId;
859         }
860     }
861   return nbvol;
862 }
863
864 /*! get the node id's of a cell.
865  * The cell is defined by it's downward connectivity id and type.
866  * @param nodeSet set of of vtk node id's to fill.
867  * @param downId downward connectivity id of the cell.
868  * @param downType type of cell.
869  */
870 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
871 {
872   _downArray[downType]->getNodeIds(downId, nodeSet);
873 }
874
875 /*! change some nodes in cell without modifying type or internal connectivity.
876  * Nodes inverse connectivity is maintained up to date.
877  * @param vtkVolId vtk id of the cell
878  * @param localClonedNodeIds map old node id to new node id.
879  */
880 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
881 {
882   vtkIdType npts = 0;
883   vtkIdType *pts; // will refer to the point id's of the face
884   this->GetCellPoints(vtkVolId, npts, pts);
885   for (int i = 0; i < npts; i++)
886     {
887       if (localClonedNodeIds.count(pts[i]))
888         {
889           vtkIdType oldpt = pts[i];
890           pts[i] = localClonedNodeIds[oldpt];
891           //MESSAGE(oldpt << " --> " << pts[i]);
892           //this->RemoveReferenceToCell(oldpt, vtkVolId);
893           //this->AddReferenceToCell(pts[i], vtkVolId);
894         }
895     }
896 }
897
898 /*! reorder the nodes of a face
899  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
900  * @param orderedNodes list of nodes to reorder (in out)
901  * @return size of the list
902  */
903 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes)
904 {
905   int vtkType = this->GetCellType(vtkVolId);
906   int cellDim = SMDS_Downward::getCellDimension(vtkType);
907   if (cellDim != 3)
908     return 0;
909   SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
910   int downVolId = this->_cellIdToDownId[vtkVolId];
911   downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
912   return orderedNodes.size();
913 }
914
915 void SMDS_UnstructuredGrid::BuildLinks()
916 {
917   // Remove the old links if they are already built
918   if (this->Links)
919     {
920     this->Links->UnRegister(this);
921     }
922
923   this->Links = SMDS_CellLinks::New();
924   this->Links->Allocate(this->GetNumberOfPoints());
925   this->Links->Register(this);
926   this->Links->BuildLinks(this, this->Connectivity);
927   this->Links->Delete();
928 }
929
930 /*! Create a volume (prism or hexahedron) by duplication of a face.
931  * Designed for use in creation of flat elements separating volume domains.
932  * A face separating two domains is shared by two volume cells.
933  * All the nodes are already created (for the two faces).
934  * Each original Node is associated to corresponding nodes in the domains.
935  * Some nodes may be duplicated for more than two domains, when domain separations intersect.
936  * In that case, even some of the nodes to use for the original face may be changed.
937  * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
938  * @param domain1: domain of the original face
939  * @param domain2: domain of the duplicated face
940  * @param originalNodes: the vtk node ids of the original face
941  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
942  * @return ok if success.
943  */
944 SMDS_MeshVolume* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
945                                                   int domain1,
946                                                   int domain2,
947                                                   std::set<int>& originalNodes,
948                                                   std::map<int, std::map<int, int> >& nodeDomains,
949                                                   std::map<int, std::map<long, int> >& nodeQuadDomains)
950 {
951   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
952   vector<vtkIdType> orderedOriginals;
953   orderedOriginals.clear();
954   set<int>::const_iterator it = originalNodes.begin();
955   for (; it != originalNodes.end(); ++it)
956     orderedOriginals.push_back(*it);
957
958   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, orderedOriginals);
959   vector<vtkIdType> orderedNodes;
960
961   switch (orderedOriginals.size())
962   {
963     case 3:
964     case 4:
965       for (int i = 0; i < nbNodes; i++)
966         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
967       for (int i = 0; i < nbNodes; i++)
968         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
969       break;
970     case 6:
971     case 8:
972       {
973         long dom1 = domain1;
974         long dom2 = domain2;
975         long dom1_2; // for nodeQuadDomains
976         if (domain1 < domain2)
977           dom1_2 = dom1 + INT_MAX * dom2;
978         else
979           dom1_2 = dom2 + INT_MAX * dom1;
980         //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
981         int ima = orderedOriginals.size();
982         int mid = orderedOriginals.size() / 2;
983         //cerr << "ima=" << ima << " mid=" << mid << endl;
984         for (int i = 0; i < mid; i++)
985           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
986         for (int i = 0; i < mid; i++)
987           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
988         for (int i = mid; i < ima; i++)
989           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
990         for (int i = mid; i < ima; i++)
991           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
992         for (int i = 0; i < mid; i++)
993           {
994             int oldId = orderedOriginals[i];
995             int newId;
996             if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
997               newId = nodeQuadDomains[oldId][dom1_2];
998             else
999               {
1000                 double *coords = this->GetPoint(oldId);
1001                 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1002                 newId = newNode->getVtkId();
1003                 std::map<long, int> emptyMap;
1004                 nodeQuadDomains[oldId] = emptyMap;
1005                 nodeQuadDomains[oldId][dom1_2] = newId;
1006               }
1007             orderedNodes.push_back(newId);
1008           }
1009       }
1010       break;
1011     default:
1012       ASSERT(0);
1013   }
1014
1015   SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1016
1017   // TODO update subshape list of elements and nodes
1018   return vol;
1019 }