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