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