Salome HOME
PR: double nodes and flat elements for ASTER calculations in progress
[modules/smesh.git] / src / SMDS / SMDS_Downward.cxx
1 /*
2  * SMDS_Downward.cxx
3  *
4  *  Created on: Jun 3, 2010
5  *      Author: prascle
6  */
7
8 #include "SMDS_Downward.hxx"
9 #include "SMDS_Mesh.hxx"
10 #include "utilities.h"
11
12 #include <vtkCellType.h>
13 #include <vtkCellLinks.h>
14
15 #include <map>
16
17 using namespace std;
18
19 // ---------------------------------------------------------------------------
20
21 vector<int> SMDS_Downward::_cellDimension;
22
23 /*! get the dimension of a cell (1,2,3 for 1D, 2D 3D) given the vtk cell type
24  *
25  * @param cellType vtk cell type @see vtkCellType.h
26  * @return 1,2 or 3
27  */
28 int SMDS_Downward::getCellDimension(unsigned char cellType)
29 {
30   if (_cellDimension.empty())
31     {
32       _cellDimension.resize(VTK_MAXTYPE + 1, 0);
33       _cellDimension[VTK_LINE] = 1;
34       _cellDimension[VTK_QUADRATIC_EDGE] = 1;
35       _cellDimension[VTK_TRIANGLE] = 2;
36       _cellDimension[VTK_QUADRATIC_TRIANGLE] = 2;
37       _cellDimension[VTK_QUAD] = 2;
38       _cellDimension[VTK_QUADRATIC_QUAD] = 2;
39       _cellDimension[VTK_TETRA] = 3;
40       _cellDimension[VTK_QUADRATIC_TETRA] = 3;
41       _cellDimension[VTK_HEXAHEDRON] = 3;
42       _cellDimension[VTK_QUADRATIC_HEXAHEDRON] = 3;
43       _cellDimension[VTK_WEDGE] = 3;
44       _cellDimension[VTK_QUADRATIC_WEDGE] = 3;
45       _cellDimension[VTK_PYRAMID] = 3;
46       _cellDimension[VTK_QUADRATIC_PYRAMID] = 3;
47     }
48   return _cellDimension[cellType];
49 }
50
51 // ---------------------------------------------------------------------------
52
53 /*! Generic constructor for all the downward connectivity structures (one per vtk cell type).
54  *  The static structure for cell dimension is set only once.
55  * @param grid unstructured grid associated to the mesh.
56  * @param nbDownCells number of downward entities associated to this vtk type of cell.
57  * @return
58  */
59 SMDS_Downward::SMDS_Downward(SMDS_UnstructuredGrid *grid, int nbDownCells) :
60   _grid(grid), _nbDownCells(nbDownCells)
61 {
62   this->_maxId = 0;
63   this->_cellIds.clear();
64   this->_cellTypes.clear();
65   if (_cellDimension.empty())
66     {
67       _cellDimension.resize(VTK_MAXTYPE + 1, 0);
68       _cellDimension[VTK_LINE] = 1;
69       _cellDimension[VTK_QUADRATIC_EDGE] = 1;
70       _cellDimension[VTK_TRIANGLE] = 2;
71       _cellDimension[VTK_QUADRATIC_TRIANGLE] = 2;
72       _cellDimension[VTK_QUAD] = 2;
73       _cellDimension[VTK_QUADRATIC_QUAD] = 2;
74       _cellDimension[VTK_TETRA] = 3;
75       _cellDimension[VTK_QUADRATIC_TETRA] = 3;
76       _cellDimension[VTK_HEXAHEDRON] = 3;
77       _cellDimension[VTK_QUADRATIC_HEXAHEDRON] = 3;
78       _cellDimension[VTK_WEDGE] = 3;
79       _cellDimension[VTK_QUADRATIC_WEDGE] = 3;
80       _cellDimension[VTK_PYRAMID] = 3;
81       _cellDimension[VTK_QUADRATIC_PYRAMID] = 3;
82     }
83 }
84
85 SMDS_Downward::~SMDS_Downward()
86 {
87 }
88
89 /*! Give or create an entry for downward connectivity structure relative to a cell.
90  * If the entry already exists, just return its id, otherwise, create it.
91  * The internal storage memory is allocated if needed.
92  * The SMDS_UnstructuredGrid::_cellIdToDownId vector is completed for vtkUnstructuredGrid cells.
93  * @param vtkId for a vtkUnstructuredGrid cell  or -1 (default) for a created downward cell.
94  * @return the rank in downward[vtkType] structure.
95  */
96 int SMDS_Downward::addCell(int vtkId)
97 {
98   int localId = -1;
99   if (vtkId >= 0)
100     localId = _grid->CellIdToDownId(vtkId);
101   if (localId >= 0)
102     return localId;
103
104   localId = this->_maxId;
105   this->_maxId++;
106   this->allocate(_maxId);
107   if (vtkId >= 0)
108     {
109       this->_vtkCellIds[localId] = vtkId;
110       _grid->setCellIdToDownId(vtkId, localId);
111     }
112   this->initCell(localId);
113   return localId;
114 }
115
116 /*! generic method do nothing. see derived methods
117  *
118  * @param cellId
119  */
120 void SMDS_Downward::initCell(int cellId)
121 {
122 }
123
124 /*! Get the number of downward entities associated to a cell (always the same for a given vtk type of cell)
125  *
126  * @param cellId not used here.
127  * @return
128  */
129 int SMDS_Downward::getNumberOfDownCells(int cellId)
130 {
131   return _nbDownCells;
132 }
133
134 /*! get a pointer on the downward entities id's associated to a cell.
135  * @see SMDS_Downward::getNumberOfDownCells for the number of downward entities.
136  * @see SMDS_Downward::getDownTypes for the vtk cell types associated to the downward entities.
137  * @param cellId index of the cell in the downward structure relative to a given vtk cell type.
138  * @return table of downward entities id's.
139  */
140 const int* SMDS_Downward::getDownCells(int cellId)
141 {
142   //ASSERT((cellId >=0) && (cellId < _maxId));
143   return &_cellIds[_nbDownCells * cellId];
144 }
145
146 /*! get a list of vtk cell types associated to downward entities of a given cell, in the same order
147  * than the downward entities id's list (@see SMDS_Downward::getDownCells).
148  *
149  * @param cellId index of the cell in the downward structure relative to a vtk cell type.
150  * @return table of downward entities types.
151  */
152 const unsigned char* SMDS_Downward::getDownTypes(int cellId)
153 {
154   return &_cellTypes[0];
155 }
156
157 /*! add a downward entity of dimension n-1 (cell or node) to a given cell.
158  * Actual implementation is done in derived methods.
159  * @param cellId index of the parent cell (dimension n) in the downward structure relative to a vtk cell type.
160  * @param lowCellId index of the children cell to add (dimension n-1)
161  * @param aType vtk cell type of the cell to add (needed to find the SMDS_Downward structure containing the cell to add).
162  */
163 void SMDS_Downward::addDownCell(int cellId, int lowCellId, unsigned char aType)
164 {
165   ASSERT(0); // must be re-implemented in derived class
166 }
167
168 /*! add a downward entity of dimension n+1 to a given cell.
169  * Actual implementation is done in derived methods.
170  * @param cellId index of the children cell (dimension n) in the downward structure relative to a vtk cell type.
171  * @param upCellId index of the parent cell to add (dimension n+1)
172  * @param aType vtk cell type of the cell to add (needed to find the SMDS_Downward structure containing the cell to add).
173  */
174 void SMDS_Downward::addUpCell(int cellId, int upCellId, unsigned char aType)
175 {
176   ASSERT(0); // must be re-implemented in derived class
177 }
178
179 int SMDS_Downward::getNodeSet(int cellId, int* nodeSet)
180 {
181   return 0;
182 }
183
184 // ---------------------------------------------------------------------------
185
186 SMDS_Down1D::SMDS_Down1D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
187   SMDS_Downward(grid, nbDownCells)
188 {
189   _upCellIdsVector.clear();
190   _upCellTypesVector.clear();
191   _upCellIds.clear();
192   _upCellTypes.clear();
193   _upCellIndex.clear();
194 }
195
196 SMDS_Down1D::~SMDS_Down1D()
197 {
198 }
199
200 /*! clear vectors used to reference 2D cells containing the edge
201  *
202  * @param cellId
203  */
204 void SMDS_Down1D::initCell(int cellId)
205 {
206   _upCellIdsVector[cellId].clear();
207   _upCellTypesVector[cellId].clear();
208 }
209
210 /*! Resize the downward connectivity storage vector if needed.
211  *
212  * @param nbElems total number of elements of the same type required
213  */
214 void SMDS_Down1D::allocate(int nbElems)
215 {
216   if (nbElems >= _vtkCellIds.size())
217     {
218       _vtkCellIds.resize(nbElems + SMDS_Mesh::chunkSize, -1);
219       _cellIds.resize(_nbDownCells * (nbElems + SMDS_Mesh::chunkSize), -1);
220       _upCellIdsVector.resize(nbElems + SMDS_Mesh::chunkSize);
221       _upCellTypesVector.resize(nbElems + SMDS_Mesh::chunkSize);
222     }
223 }
224
225 void SMDS_Down1D::compactStorage()
226 {
227   _cellIds.resize(_nbDownCells * _maxId);
228   _vtkCellIds.resize(_maxId);
229
230   int sizeUpCells = 0;
231   for (int i = 0; i < _maxId; i++)
232     sizeUpCells += _upCellIdsVector[i].size();
233   _upCellIds.resize(sizeUpCells, -1);
234   _upCellTypes.resize(sizeUpCells);
235   _upCellIndex.resize(_maxId + 1, -1); // id and types of rank i correspond to [ _upCellIndex[i], _upCellIndex[i+1] [
236   int current = 0;
237   for (int i = 0; i < _maxId; i++)
238     {
239       _upCellIndex[i] = current;
240       for (int j = 0; j < _upCellIdsVector[i].size(); j++)
241         {
242           _upCellIds[current] = _upCellIdsVector[i][j];
243           _upCellTypes[current] = _upCellTypesVector[i][j];
244           current++;
245         }
246     }
247   _upCellIndex[_maxId] = current;
248
249   _upCellIdsVector.clear();
250   _upCellTypesVector.clear();
251 }
252
253 void SMDS_Down1D::addUpCell(int cellId, int upCellId, unsigned char aType)
254 {
255   //ASSERT((cellId >=0) && (cellId < _maxId));
256   int nbFaces = _upCellIdsVector[cellId].size();
257   for (int i = 0; i < nbFaces; i++)
258     {
259       if ((_upCellIdsVector[cellId][i] == upCellId) && (_upCellTypesVector[cellId][i] == aType))
260         {
261           return; // already done
262         }
263     }
264   _upCellIdsVector[cellId].push_back(upCellId);
265   _upCellTypesVector[cellId].push_back(aType);
266 }
267
268 int SMDS_Down1D::getNumberOfUpCells(int cellId)
269 {
270   //ASSERT((cellId >=0) && (cellId < _maxId));
271   return _upCellIndex[cellId + 1] - _upCellIndex[cellId];
272 }
273
274 const int* SMDS_Down1D::getUpCells(int cellId)
275 {
276   //ASSERT((cellId >=0) && (cellId < _maxId));
277   return &_upCellIds[_upCellIndex[cellId]];
278 }
279
280 const unsigned char* SMDS_Down1D::getUpTypes(int cellId)
281 {
282   //ASSERT((cellId >=0) && (cellId < _maxId));
283   return &_upCellTypes[_upCellIndex[cellId]];
284 }
285
286 void SMDS_Down1D::getNodeIds(int cellId, std::set<int>& nodeSet)
287 {
288   for (int i = 0; i < _nbDownCells; i++)
289     nodeSet.insert(_cellIds[_nbDownCells * cellId + i]);
290 }
291
292 int SMDS_Down1D::getNodeSet(int cellId, int* nodeSet)
293 {
294   for (int i = 0; i < _nbDownCells; i++)
295     nodeSet[i] = _cellIds[_nbDownCells * cellId + i];
296   return _nbDownCells;
297 }
298
299 void SMDS_Down1D::setNodes(int cellId, int vtkId)
300 {
301   vtkIdType npts = 0;
302   vtkIdType *pts; // will refer to the point id's of the face
303   _grid->GetCellPoints(vtkId, npts, pts);
304   // MESSAGE(vtkId << " " << npts << "  " << _nbDownCells);
305   //ASSERT(npts == _nbDownCells);
306   for (int i = 0; i < npts; i++)
307     {
308       _cellIds[_nbDownCells * cellId + i] = pts[i];
309     }
310 }
311
312 void SMDS_Down1D::setNodes(int cellId, const int* nodeIds)
313 {
314   //ASSERT(nodeIds.size() == _nbDownCells);
315   for (int i = 0; i < _nbDownCells; i++)
316     {
317       _cellIds[_nbDownCells * cellId + i] = nodeIds[i];
318     }
319 }
320
321 /*! Build the list of vtkUnstructuredGrid cells containing the edge.
322  * We keep in the list the cells that contains all the nodes, we keep only volumes and faces.
323  * @param cellId id of the edge in the downward structure
324  * @param vtkIds vector of vtk id's
325  * @return number of vtk cells (size of vector)
326  */
327 int SMDS_Down1D::computeVtkCells(int cellId, std::vector<int>& vtkIds)
328 {
329   vtkIds.clear();
330
331   // --- find all the cells the points belong to, and how many of the points belong to a given cell
332
333   int *pts = &_cellIds[_nbDownCells * cellId];
334   int ncells = this->computeVtkCells(pts, vtkIds);
335   return ncells;
336 }
337
338 /*! Build the list of vtkUnstructuredGrid cells containing the edge.
339  *
340  * @param pts list of points id's defining an edge
341  * @param vtkIds vector of vtk id's
342  * @return number of vtk cells (size of vector)
343  */
344 int SMDS_Down1D::computeVtkCells(int *pts, std::vector<int>& vtkIds)
345 {
346
347   // --- find all the cells the points belong to, and how many of the points belong to a given cell
348
349   int cellIds[1000];
350   int cellCnt[1000];
351   int cnt = 0;
352   for (int i = 0; i < _nbDownCells; i++)
353     {
354       vtkIdType point = pts[i];
355       int numCells = _grid->GetLinks()->GetNcells(point);
356       vtkIdType *cells = _grid->GetLinks()->GetCells(point);
357       for (int j = 0; j < numCells; j++)
358         {
359           int vtkCellId = cells[j];
360           bool found = false;
361           for (int k = 0; k < cnt; k++)
362             {
363               if (cellIds[k] == vtkCellId)
364                 {
365                   cellCnt[k] += 1;
366                   found = true;
367                   break;
368                 }
369             }
370           if (!found)
371             {
372               cellIds[cnt] = vtkCellId;
373               cellCnt[cnt] = 1;
374               // TODO ASSERT(cnt<1000);
375               cnt++;
376             }
377         }
378     }
379
380   // --- find the face and volume cells: they contains all the points and are of type volume or face
381
382   int ncells = 0;
383   for (int i = 0; i < cnt; i++)
384     {
385       if (cellCnt[i] == _nbDownCells)
386         {
387           int vtkElemId = cellIds[i];
388           int vtkType = _grid->GetCellType(vtkElemId);
389           if (SMDS_Downward::getCellDimension(vtkType) > 1)
390             {
391               vtkIds.push_back(vtkElemId);
392               ncells++;
393             }
394         }
395     }
396
397   return ncells;
398 }
399
400 /*! Build the list of downward faces from a list of vtk cells.
401  *
402  * @param cellId id of the edge in the downward structure
403  * @param vtkIds vector of vtk id's
404  * @param downFaces vector of face id's in downward structures
405  * @param downTypes vector of face types
406  * @return number of downward faces
407  */
408 int SMDS_Down1D::computeFaces(int cellId, int* vtkIds, int nbcells, int* downFaces, unsigned char* downTypes)
409 {
410   int *pts = &_cellIds[_nbDownCells * cellId];
411   int nbFaces = this->computeFaces(pts, vtkIds, nbcells, downFaces, downTypes);
412   return nbFaces;
413 }
414
415 /*! Build the list of downward faces from a list of vtk cells.
416  *
417  * @param pts list of points id's defining an edge
418  * @param vtkIds vector of vtk id's
419  * @param downFaces vector of face id's in downward structures
420  * @param downTypes vector of face types
421  * @return number of downward faces
422  */
423 int SMDS_Down1D::computeFaces(int* pts, int* vtkIds, int nbcells, int* downFaces, unsigned char* downTypes)
424 {
425   int cnt = 0;
426   for (int i = 0; i < nbcells; i++)
427     {
428       int vtkId = vtkIds[i];
429       int vtkType = _grid->GetCellType(vtkId);
430       if (SMDS_Downward::getCellDimension(vtkType) == 2)
431         {
432           int faceId = _grid->CellIdToDownId(vtkId);
433           downFaces[cnt] = faceId;
434           downTypes[cnt] = vtkType;
435           cnt++;
436         }
437       else if (SMDS_Downward::getCellDimension(vtkType) == 3)
438         {
439           int volId = _grid->CellIdToDownId(vtkId);
440           SMDS_Downward * downvol = _grid->getDownArray(vtkType);
441           //const int *downIds = downvol->getDownCells(volId);
442           const unsigned char* downTypesVol = downvol->getDownTypes(volId);
443           int nbFaces = downvol->getNumberOfDownCells(volId);
444           const int* faceIds = downvol->getDownCells(volId);
445           for (int n = 0; n < nbFaces; n++)
446             {
447               SMDS_Down2D *downFace = static_cast<SMDS_Down2D*> (_grid->getDownArray(downTypesVol[n]));
448               bool isInFace = downFace->isInFace(faceIds[n], pts, _nbDownCells);
449               if (isInFace)
450                 {
451                   bool alreadySet = false;
452                   for (int k = 0; k < cnt; k++)
453                     if (faceIds[n] == downFaces[k])
454                       {
455                         alreadySet = true;
456                         break;
457                       }
458                   if (!alreadySet)
459                     {
460                       downFaces[cnt] = faceIds[n];
461                       downTypes[cnt] = downTypesVol[n];
462                       cnt++;
463                     }
464                 }
465             }
466         }
467     }
468   return cnt;
469 }
470
471 // ---------------------------------------------------------------------------
472
473 SMDS_Down2D::SMDS_Down2D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
474   SMDS_Downward(grid, nbDownCells)
475 {
476   _upCellIds.clear();
477   _upCellTypes.clear();
478   _tempNodes.clear();
479   _nbNodes = 0;
480 }
481
482 SMDS_Down2D::~SMDS_Down2D()
483 {
484 }
485
486 int SMDS_Down2D::getNumberOfUpCells(int cellId)
487 {
488   int nbup = 0;
489   if (_upCellIds[2 * cellId] >= 0)
490     nbup++;
491   if (_upCellIds[2 * cellId + 1] >= 0)
492     nbup++;
493   return nbup;
494 }
495
496 const int* SMDS_Down2D::getUpCells(int cellId)
497 {
498   //ASSERT((cellId >=0) && (cellId < _maxId));
499   return &_upCellIds[2 * cellId];
500 }
501
502 const unsigned char* SMDS_Down2D::getUpTypes(int cellId)
503 {
504   //ASSERT((cellId >=0) && (cellId < _maxId));
505   return &_upCellTypes[2 * cellId];
506 }
507
508 void SMDS_Down2D::getNodeIds(int cellId, std::set<int>& nodeSet)
509 {
510   for (int i = 0; i < _nbDownCells; i++)
511     {
512       int downCellId = _cellIds[_nbDownCells * cellId + i];
513       unsigned char cellType = _cellTypes[i];
514       this->_grid->getDownArray(cellType)->getNodeIds(downCellId, nodeSet);
515     }
516 }
517
518 /*! Find in vtkUnstructuredGrid the volumes containing a face already stored in vtkUnstructuredGrid.
519  * Search the volumes containing a face, to store the info in SMDS_Down2D for later uses
520  * with SMDS_Down2D::getUpCells and SMDS_Down2D::getUpTypes.
521  * A face belongs to 0, 1 or 2 volumes, identified by their id in vtkUnstructuredGrid.
522  * @param cellId the face cell id in vkUnstructuredGrid
523  * @param ids a couple of vtkId, initialized at -1 (no parent volume)
524  * @return number of volumes (0, 1 or 2)
525  */
526 int SMDS_Down2D::computeVolumeIds(int cellId, int* ids)
527 {
528   // --- find point id's of the face
529
530   vtkIdType npts = 0;
531   vtkIdType *pts; // will refer to the point id's of the face
532   _grid->GetCellPoints(cellId, npts, pts);
533   vector<int> nodes;
534   for (int i = 0; i < npts; i++)
535     nodes.push_back(pts[i]);
536   int nvol = this->computeVolumeIdsFromNodesFace(&nodes[0], npts, ids);
537   return nvol;
538 }
539
540 /*! Find in vtkUnstructuredGrid the volumes containing a face described by it's nodes
541  * Search the volumes containing a face, to store the info in SMDS_Down2D for later uses
542  * with SMDS_Down2D::getUpCells and SMDS_Down2D::getUpTypes.
543  * A face belongs to 0, 1 or 2 volumes, identified by their id in vtkUnstructuredGrid.
544  * @param faceByNodes
545  * @param ids a couple of vtkId, initialized at -1 (no parent volume)
546  * @return number of volumes (0, 1 or 2)
547  */
548 int SMDS_Down2D::computeVolumeIds(ElemByNodesType& faceByNodes, int* ids)
549 {
550   int nvol = this->computeVolumeIdsFromNodesFace(&faceByNodes.nodeIds[0], faceByNodes.nbNodes, ids);
551   return nvol;
552 }
553
554 /*! Find in vtkUnstructuredGrid the volumes containing a face described by it's nodes
555  * Search the volumes containing a face, to store the info in SMDS_Down2D for later uses
556  * with SMDS_Down2D::getUpCells and SMDS_Down2D::getUpTypes.
557  * A face belongs to 0, 1 or 2 volumes, identified by their id in vtkUnstructuredGrid.
558  * @param pts array of vtk node id's
559  * @param npts number of nodes
560  * @param ids
561  * @return number of volumes (0, 1 or 2)
562  */
563 int SMDS_Down2D::computeVolumeIdsFromNodesFace(int* pts, int npts, int* ids)
564 {
565
566   // --- find all the cells the points belong to, and how many of the points belong to a given cell
567
568   int cellIds[1000];
569   int cellCnt[1000];
570   int cnt = 0;
571   for (int i = 0; i < npts; i++)
572     {
573       vtkIdType point = pts[i];
574       int numCells = _grid->GetLinks()->GetNcells(point);
575       //MESSAGE("cells pour " << i << " " << numCells);
576       vtkIdType *cells = _grid->GetLinks()->GetCells(point);
577       for (int j = 0; j < numCells; j++)
578         {
579           int vtkCellId = cells[j];
580           bool found = false;
581           for (int k = 0; k < cnt; k++)
582             {
583               if (cellIds[k] == vtkCellId)
584                 {
585                   cellCnt[k] += 1;
586                   found = true;
587                   break;
588                 }
589             }
590           if (!found)
591             {
592               cellIds[cnt] = vtkCellId;
593               cellCnt[cnt] = 1;
594               // TODO ASSERT(cnt<1000);
595               cnt++;
596             }
597         }
598     }
599
600   // --- find the volume cells: they contains all the points and are of type volume
601
602   int nvol = 0;
603   for (int i = 0; i < cnt; i++)
604     {
605       //MESSAGE("cell " << cellIds[i] << " points " << cellCnt[i]);
606       if (cellCnt[i] == npts)
607         {
608           int vtkElemId = cellIds[i];
609           int vtkType = _grid->GetCellType(vtkElemId);
610           if (SMDS_Downward::getCellDimension(vtkType) == 3)
611             {
612               ids[nvol] = vtkElemId; // store the volume id in given vector
613               nvol++;
614             }
615         }
616       if (nvol == 2)
617         break;
618     }
619
620   return nvol;
621 }
622
623 void SMDS_Down2D::setTempNodes(int cellId, int vtkId)
624 {
625   vtkIdType npts = 0;
626   vtkIdType *pts; // will refer to the point id's of the face
627   _grid->GetCellPoints(vtkId, npts, pts);
628   // MESSAGE(vtkId << " " << npts << "  " << _nbNodes);
629   //ASSERT(npts == _nbNodes);
630   for (int i = 0; i < npts; i++)
631     {
632       _tempNodes[_nbNodes * cellId + i] = pts[i];
633     }
634 }
635
636 void SMDS_Down2D::setTempNodes(int cellId, ElemByNodesType& faceByNodes)
637 {
638   for (int i = 0; i < faceByNodes.nbNodes; i++)
639     _tempNodes[_nbNodes * cellId + i] = faceByNodes.nodeIds[i];
640 }
641
642 /*! Find if all the nodes belongs to the face.
643  *
644  * @param cellId the face cell Id
645  * @param nodeSet set of node id's to be found in the face list of nodes
646  * @return
647  */
648 bool SMDS_Down2D::isInFace(int cellId, int *pts, int npts)
649 {
650   int nbFound = 0;
651   int *nodes = &_tempNodes[_nbNodes * cellId];
652   for (int j = 0; j < npts; j++)
653     {
654       int point = pts[j];
655       for (int i = 0; i < _nbNodes; i++)
656         {
657           if (nodes[i] == point)
658             {
659               nbFound++;
660               break;
661             }
662         }
663     }
664   return (nbFound == npts);
665 }
666
667 /*! Resize the downward connectivity storage vector if needed.
668  *
669  * @param nbElems total number of elements of the same type required
670  */
671 void SMDS_Down2D::allocate(int nbElems)
672 {
673   if (nbElems >= _vtkCellIds.size())
674     {
675       _cellIds.resize(_nbDownCells * (nbElems + SMDS_Mesh::chunkSize), -1);
676       _vtkCellIds.resize(nbElems + SMDS_Mesh::chunkSize, -1);
677       _upCellIds.resize(2 * (nbElems + SMDS_Mesh::chunkSize), -1);
678       _upCellTypes.resize(2 * (nbElems + SMDS_Mesh::chunkSize), -1);
679       _tempNodes.resize(_nbNodes * (nbElems + SMDS_Mesh::chunkSize), -1);
680     }
681 }
682
683 void SMDS_Down2D::compactStorage()
684 {
685   _cellIds.resize(_nbDownCells * _maxId);
686   _upCellIds.resize(2 * _maxId);
687   _upCellTypes.resize(2 * _maxId);
688   _vtkCellIds.resize(_maxId);
689   _tempNodes.clear();
690 }
691
692 void SMDS_Down2D::addUpCell(int cellId, int upCellId, unsigned char aType)
693 {
694   //ASSERT((cellId >=0)&& (cellId < _maxId));
695   int *vols = &_upCellIds[2 * cellId];
696   unsigned char *types = &_upCellTypes[2 * cellId];
697   for (int i = 0; i < 2; i++)
698     {
699       if (vols[i] < 0)
700         {
701           vols[i] = upCellId; // use non affected volume
702           types[i] = aType;
703           return;
704         }
705       if ((vols[i] == upCellId) && (types[i] == aType)) // already done
706         return;
707     }
708   ASSERT(0);
709 }
710
711 int SMDS_Down2D::getNodeSet(int cellId, int* nodeSet)
712 {
713   for (int i = 0; i < _nbNodes; i++)
714     nodeSet[i] = _tempNodes[_nbNodes * cellId + i];
715   return _nbNodes;
716 }
717
718 int SMDS_Down2D::FindEdgeByNodes(int cellId, ElemByNodesType& edgeByNodes)
719 {
720   int *edges = &_cellIds[_nbDownCells * cellId];
721   for (int i = 0; i < _nbDownCells; i++)
722     {
723       if ((edges[i] >= 0) && (edgeByNodes.vtkType == _cellTypes[i]))
724         {
725           int nodeSet[3];
726           int npts = this->_grid->getDownArray(edgeByNodes.vtkType)->getNodeSet(edges[i], nodeSet);
727           bool found = false;
728           for (int j = 0; j < npts; j++)
729             {
730               int point = edgeByNodes.nodeIds[j];
731               found = false;
732               for (int k = 0; k < npts; k++)
733                 {
734                   if (nodeSet[k] == point)
735                     {
736                       found = true;
737                       break;
738                     }
739                 }
740               if (!found)
741                 break;
742             }
743           if (found)
744             return edges[i];
745         }
746     }
747   return -1;
748 }
749
750 // ---------------------------------------------------------------------------
751
752 SMDS_Down3D::SMDS_Down3D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
753   SMDS_Downward(grid, nbDownCells)
754 {
755 }
756
757 SMDS_Down3D::~SMDS_Down3D()
758 {
759 }
760
761 void SMDS_Down3D::allocate(int nbElems)
762 {
763   if (nbElems >= _vtkCellIds.size())
764     {
765       _cellIds.resize(_nbDownCells * (nbElems + SMDS_Mesh::chunkSize), -1);
766       _vtkCellIds.resize(nbElems + SMDS_Mesh::chunkSize, -1);
767     }
768 }
769
770 void SMDS_Down3D::compactStorage()
771 {
772   // nothing to do, size was known before
773 }
774
775 int SMDS_Down3D::getNumberOfUpCells(int cellId)
776 {
777   return 0;
778 }
779
780 const int* SMDS_Down3D::getUpCells(int cellId)
781 {
782   return 0;
783 }
784
785 const unsigned char* SMDS_Down3D::getUpTypes(int cellId)
786 {
787   return 0;
788 }
789
790 void SMDS_Down3D::getNodeIds(int cellId, std::set<int>& nodeSet)
791 {
792   int vtkId = this->_vtkCellIds[cellId];
793   vtkIdType npts = 0;
794   vtkIdType *nodes; // will refer to the point id's of the volume
795   _grid->GetCellPoints(vtkId, npts, nodes);
796   for (int i = 0; i < npts; i++)
797     nodeSet.insert(nodes[i]);
798 }
799
800 int SMDS_Down3D::FindFaceByNodes(int cellId, ElemByNodesType& faceByNodes)
801 {
802   int *faces = &_cellIds[_nbDownCells * cellId];
803   int faceNodeSet[10];
804   int npoints = 0;
805
806   for (int i = 0; i < _nbDownCells; i++)
807     {
808       if ((faces[i] >= 0) && (faceByNodes.vtkType == _cellTypes[i]))
809         {
810           if (npoints == 0)
811             {
812               for (int j = 0; j < faceByNodes.nbNodes; j++)
813                 faceNodeSet[j] = faceByNodes.nodeIds[j];
814               npoints = faceByNodes.nbNodes;
815             }
816
817           int nodeSet[10];
818           int npts = this->_grid->getDownArray(faceByNodes.vtkType)->getNodeSet(faces[i], nodeSet);
819           if (npts != npoints)
820             continue; // skip this face
821           bool found = false;
822           for (int j = 0; j < npts; j++)
823             {
824               int point = faceByNodes.nodeIds[j];
825               found = false;
826               for (int k = 0; k < npts; k++)
827                 {
828                   if (nodeSet[k] == point)
829                     {
830                       found = true;
831                       break; // point j is in the 2 faces, skip remaining k values
832                     }
833                 }
834               if (!found)
835                 break; // point j is not in the 2 faces, skip the remaining tests
836             }
837           if (found)
838             return faces[i];
839         }
840     }
841   return -1;
842 }
843
844 // ---------------------------------------------------------------------------
845
846 SMDS_DownEdge::SMDS_DownEdge(SMDS_UnstructuredGrid *grid) :
847   SMDS_Down1D(grid, 2)
848 {
849   _cellTypes.push_back(VTK_VERTEX);
850   _cellTypes.push_back(VTK_VERTEX);
851 }
852
853 SMDS_DownEdge::~SMDS_DownEdge()
854 {
855 }
856
857 // ---------------------------------------------------------------------------
858
859 SMDS_DownQuadEdge::SMDS_DownQuadEdge(SMDS_UnstructuredGrid *grid) :
860   SMDS_Down1D(grid, 3)
861 {
862   _cellTypes.push_back(VTK_VERTEX);
863   _cellTypes.push_back(VTK_VERTEX);
864   _cellTypes.push_back(VTK_VERTEX);
865 }
866
867 SMDS_DownQuadEdge::~SMDS_DownQuadEdge()
868 {
869 }
870
871 // ---------------------------------------------------------------------------
872
873 SMDS_DownTriangle::SMDS_DownTriangle(SMDS_UnstructuredGrid *grid) :
874   SMDS_Down2D(grid, 3)
875 {
876   _cellTypes.push_back(VTK_LINE);
877   _cellTypes.push_back(VTK_LINE);
878   _cellTypes.push_back(VTK_LINE);
879   _nbNodes = 3;
880 }
881
882 SMDS_DownTriangle::~SMDS_DownTriangle()
883 {
884 }
885
886 void SMDS_DownTriangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
887 {
888   int *nodes = &_tempNodes[_nbNodes * cellId];
889   edgesWithNodes.nbElems = 3;
890
891   edgesWithNodes.elems[0].nodeIds[0] = nodes[0];
892   edgesWithNodes.elems[0].nodeIds[1] = nodes[1];
893   edgesWithNodes.elems[0].nbNodes = 2;
894   edgesWithNodes.elems[0].vtkType = VTK_LINE;
895
896   edgesWithNodes.elems[1].nodeIds[0] = nodes[1];
897   edgesWithNodes.elems[1].nodeIds[1] = nodes[2];
898   edgesWithNodes.elems[1].nbNodes = 2;
899   edgesWithNodes.elems[1].vtkType = VTK_LINE;
900
901   edgesWithNodes.elems[2].nodeIds[0] = nodes[2];
902   edgesWithNodes.elems[2].nodeIds[1] = nodes[0];
903   edgesWithNodes.elems[2].nbNodes = 2;
904   edgesWithNodes.elems[2].vtkType = VTK_LINE;
905 }
906
907 void SMDS_DownTriangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
908 {
909   //ASSERT((cellId >=0)&& (cellId < _maxId));
910   //ASSERT(aType == VTK_LINE);
911   int *faces = &_cellIds[_nbDownCells * cellId];
912   for (int i = 0; i < _nbDownCells; i++)
913     {
914       if (faces[i] < 0)
915         {
916           faces[i] = lowCellId;
917           return;
918         }
919       if (faces[i] == lowCellId)
920         return;
921     }
922   ASSERT(0);
923 }
924
925 // ---------------------------------------------------------------------------
926
927 SMDS_DownQuadTriangle::SMDS_DownQuadTriangle(SMDS_UnstructuredGrid *grid) :
928   SMDS_Down2D(grid, 3)
929 {
930   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
931   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
932   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
933   _nbNodes = 6;
934 }
935
936 SMDS_DownQuadTriangle::~SMDS_DownQuadTriangle()
937 {
938 }
939
940 void SMDS_DownQuadTriangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
941 {
942   int *nodes = &_tempNodes[_nbNodes * cellId];
943   edgesWithNodes.nbElems = 3;
944
945   edgesWithNodes.elems[0].nodeIds[0] = nodes[0];
946   edgesWithNodes.elems[0].nodeIds[1] = nodes[1];
947   edgesWithNodes.elems[0].nodeIds[2] = nodes[3];
948   edgesWithNodes.elems[0].nbNodes = 3;
949   edgesWithNodes.elems[0].vtkType = VTK_QUADRATIC_EDGE;
950
951   edgesWithNodes.elems[1].nodeIds[0] = nodes[1];
952   edgesWithNodes.elems[1].nodeIds[1] = nodes[2];
953   edgesWithNodes.elems[1].nodeIds[2] = nodes[4];
954   edgesWithNodes.elems[1].nbNodes = 3;
955   edgesWithNodes.elems[1].vtkType = VTK_QUADRATIC_EDGE;
956
957   edgesWithNodes.elems[2].nodeIds[0] = nodes[2];
958   edgesWithNodes.elems[2].nodeIds[1] = nodes[0];
959   edgesWithNodes.elems[2].nodeIds[2] = nodes[5];
960   edgesWithNodes.elems[2].nbNodes = 3;
961   edgesWithNodes.elems[2].vtkType = VTK_QUADRATIC_EDGE;
962 }
963
964 void SMDS_DownQuadTriangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
965 {
966   //ASSERT((cellId >=0)&& (cellId < _maxId));
967   //ASSERT(aType == VTK_QUADRATIC_EDGE);
968   int *edges = &_cellIds[_nbDownCells * cellId];
969   for (int i = 0; i < _nbDownCells; i++)
970     {
971       if (edges[i] < 0)
972         {
973           edges[i] = lowCellId;
974           return;
975         }
976       if (edges[i] == lowCellId)
977         return;
978     }
979   ASSERT(0);
980 }
981
982 // ---------------------------------------------------------------------------
983
984 SMDS_DownQuadrangle::SMDS_DownQuadrangle(SMDS_UnstructuredGrid *grid) :
985   SMDS_Down2D(grid, 4)
986 {
987   _cellTypes.push_back(VTK_LINE);
988   _cellTypes.push_back(VTK_LINE);
989   _cellTypes.push_back(VTK_LINE);
990   _cellTypes.push_back(VTK_LINE);
991   _nbNodes = 4;
992 }
993
994 SMDS_DownQuadrangle::~SMDS_DownQuadrangle()
995 {
996 }
997
998 void SMDS_DownQuadrangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
999 {
1000   int *nodes = &_tempNodes[_nbNodes * cellId];
1001   edgesWithNodes.nbElems = 4;
1002
1003   edgesWithNodes.elems[0].nodeIds[0] = nodes[0];
1004   edgesWithNodes.elems[0].nodeIds[1] = nodes[1];
1005   edgesWithNodes.elems[0].nbNodes = 2;
1006   edgesWithNodes.elems[0].vtkType = VTK_LINE;
1007
1008   edgesWithNodes.elems[1].nodeIds[0] = nodes[1];
1009   edgesWithNodes.elems[1].nodeIds[1] = nodes[2];
1010   edgesWithNodes.elems[1].nbNodes = 2;
1011   edgesWithNodes.elems[1].vtkType = VTK_LINE;
1012
1013   edgesWithNodes.elems[2].nodeIds[0] = nodes[2];
1014   edgesWithNodes.elems[2].nodeIds[1] = nodes[3];
1015   edgesWithNodes.elems[2].nbNodes = 2;
1016   edgesWithNodes.elems[2].vtkType = VTK_LINE;
1017
1018   edgesWithNodes.elems[3].nodeIds[0] = nodes[3];
1019   edgesWithNodes.elems[3].nodeIds[1] = nodes[0];
1020   edgesWithNodes.elems[3].nbNodes = 2;
1021   edgesWithNodes.elems[3].vtkType = VTK_LINE;
1022 }
1023
1024 void SMDS_DownQuadrangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
1025 {
1026   //ASSERT((cellId >=0)&& (cellId < _maxId));
1027   //ASSERT(aType == VTK_LINE);
1028   int *faces = &_cellIds[_nbDownCells * cellId];
1029   for (int i = 0; i < _nbDownCells; i++)
1030     {
1031       if (faces[i] < 0)
1032         {
1033           faces[i] = lowCellId;
1034           return;
1035         }
1036       if (faces[i] == lowCellId)
1037         return;
1038     }
1039   ASSERT(0);
1040 }
1041
1042 // ---------------------------------------------------------------------------
1043
1044 SMDS_DownQuadQuadrangle::SMDS_DownQuadQuadrangle(SMDS_UnstructuredGrid *grid) :
1045   SMDS_Down2D(grid, 4)
1046 {
1047   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
1048   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
1049   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
1050   _cellTypes.push_back(VTK_QUADRATIC_EDGE);
1051   _nbNodes = 8;
1052 }
1053
1054 SMDS_DownQuadQuadrangle::~SMDS_DownQuadQuadrangle()
1055 {
1056 }
1057
1058 void SMDS_DownQuadQuadrangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
1059 {
1060   int *nodes = &_tempNodes[_nbNodes * cellId];
1061   edgesWithNodes.nbElems = 4;
1062
1063   edgesWithNodes.elems[0].nodeIds[0] = nodes[0];
1064   edgesWithNodes.elems[0].nodeIds[1] = nodes[1];
1065   edgesWithNodes.elems[0].nodeIds[2] = nodes[4];
1066   edgesWithNodes.elems[0].nbNodes = 3;
1067   edgesWithNodes.elems[0].vtkType = VTK_QUADRATIC_EDGE;
1068
1069   edgesWithNodes.elems[1].nodeIds[0] = nodes[1];
1070   edgesWithNodes.elems[1].nodeIds[1] = nodes[2];
1071   edgesWithNodes.elems[1].nodeIds[2] = nodes[5];
1072   edgesWithNodes.elems[1].nbNodes = 3;
1073   edgesWithNodes.elems[1].vtkType = VTK_QUADRATIC_EDGE;
1074
1075   edgesWithNodes.elems[2].nodeIds[0] = nodes[2];
1076   edgesWithNodes.elems[2].nodeIds[1] = nodes[3];
1077   edgesWithNodes.elems[2].nodeIds[2] = nodes[6];
1078   edgesWithNodes.elems[2].nbNodes = 3;
1079   edgesWithNodes.elems[2].vtkType = VTK_QUADRATIC_EDGE;
1080
1081   edgesWithNodes.elems[3].nodeIds[0] = nodes[3];
1082   edgesWithNodes.elems[3].nodeIds[1] = nodes[0];
1083   edgesWithNodes.elems[3].nodeIds[2] = nodes[7];
1084   edgesWithNodes.elems[3].nbNodes = 3;
1085   edgesWithNodes.elems[3].vtkType = VTK_QUADRATIC_EDGE;
1086 }
1087
1088 void SMDS_DownQuadQuadrangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
1089 {
1090   //ASSERT((cellId >=0)&& (cellId < _maxId));
1091   //ASSERT(aType == VTK_QUADRATIC_EDGE);
1092   int *faces = &_cellIds[_nbDownCells * cellId];
1093   for (int i = 0; i < _nbDownCells; i++)
1094     {
1095       if (faces[i] < 0)
1096         {
1097           faces[i] = lowCellId;
1098           return;
1099         }
1100       if (faces[i] == lowCellId)
1101         return;
1102     }
1103   ASSERT(0);
1104 }
1105
1106 // ---------------------------------------------------------------------------
1107
1108 SMDS_DownTetra::SMDS_DownTetra(SMDS_UnstructuredGrid *grid) :
1109   SMDS_Down3D(grid, 4)
1110 {
1111   _cellTypes.push_back(VTK_TRIANGLE);
1112   _cellTypes.push_back(VTK_TRIANGLE);
1113   _cellTypes.push_back(VTK_TRIANGLE);
1114   _cellTypes.push_back(VTK_TRIANGLE);
1115 }
1116
1117 SMDS_DownTetra::~SMDS_DownTetra()
1118 {
1119 }
1120
1121 void SMDS_DownTetra::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1122 {
1123   set<int> setNodes;
1124   setNodes.clear();
1125   for (int i = 0; i < orderedNodes.size(); i++)
1126     setNodes.insert(orderedNodes[i]);
1127   //MESSAGE("cellId = " << cellId);
1128
1129   vtkIdType npts = 0;
1130   vtkIdType *nodes; // will refer to the point id's of the volume
1131   _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
1132
1133   set<int> tofind;
1134   int ids[12] = { 0, 1, 2,  0, 3, 1,  2, 3, 0,   1, 3, 2 };
1135 //int ids[12] = { 2, 1, 0,  1, 3, 0,  0, 3, 2,   2, 3, 1 };
1136   for (int k = 0; k < 4; k++)
1137     {
1138       tofind.clear();
1139       for (int i = 0; i < 3; i++)
1140         tofind.insert(nodes[ids[3 * k + i]]);
1141       if (setNodes == tofind)
1142         {
1143           for (int i = 0; i < 3; i++)
1144             orderedNodes[i] = nodes[ids[3 * k + i]];
1145           return;
1146         }
1147     }
1148   MESSAGE("=== Problem volume " << _vtkCellIds[cellId] << " " << _grid->_mesh->fromVtkToSmds(_vtkCellIds[cellId]));
1149   MESSAGE(orderedNodes[0] << " " << orderedNodes[1] << " " << orderedNodes[2]);
1150   MESSAGE(nodes[0] << " " << nodes[1] << " " << nodes[2] << " " << nodes[3]);
1151 }
1152
1153 void SMDS_DownTetra::addDownCell(int cellId, int lowCellId, unsigned char aType)
1154 {
1155   //ASSERT((cellId >=0)&& (cellId < _maxId));
1156   //ASSERT(aType == VTK_TRIANGLE);
1157   int *faces = &_cellIds[_nbDownCells * cellId];
1158   for (int i = 0; i < _nbDownCells; i++)
1159     {
1160       if (faces[i] < 0)
1161         {
1162           faces[i] = lowCellId;
1163           return;
1164         }
1165       if (faces[i] == lowCellId)
1166         return;
1167     }
1168   ASSERT(0);
1169 }
1170
1171 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1172  * The linear tetrahedron is defined by four points.
1173  * @see vtkTetra.h in Filtering.
1174  * @param cellId volumeId in vtkUnstructuredGrid
1175  * @param facesWithNodes vector of face descriptors to be filled
1176  */
1177 void SMDS_DownTetra::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1178 {
1179   // --- find point id's of the volume
1180
1181   vtkIdType npts = 0;
1182   vtkIdType *nodes; // will refer to the point id's of the volume
1183   _grid->GetCellPoints(cellId, npts, nodes);
1184
1185   // --- create all the ordered list of node id's for each face
1186
1187   facesWithNodes.nbElems = 4;
1188
1189   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1190   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
1191   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
1192   facesWithNodes.elems[0].nbNodes = 3;
1193   facesWithNodes.elems[0].vtkType = VTK_TRIANGLE;
1194
1195   facesWithNodes.elems[1].nodeIds[0] = nodes[0];
1196   facesWithNodes.elems[1].nodeIds[1] = nodes[1];
1197   facesWithNodes.elems[1].nodeIds[2] = nodes[3];
1198   facesWithNodes.elems[1].nbNodes = 3;
1199   facesWithNodes.elems[1].vtkType = VTK_TRIANGLE;
1200
1201   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
1202   facesWithNodes.elems[2].nodeIds[1] = nodes[2];
1203   facesWithNodes.elems[2].nodeIds[2] = nodes[3];
1204   facesWithNodes.elems[2].nbNodes = 3;
1205   facesWithNodes.elems[2].vtkType = VTK_TRIANGLE;
1206
1207   facesWithNodes.elems[3].nodeIds[0] = nodes[1];
1208   facesWithNodes.elems[3].nodeIds[1] = nodes[2];
1209   facesWithNodes.elems[3].nodeIds[2] = nodes[3];
1210   facesWithNodes.elems[3].nbNodes = 3;
1211   facesWithNodes.elems[3].vtkType = VTK_TRIANGLE;
1212 }
1213
1214 // ---------------------------------------------------------------------------
1215
1216 SMDS_DownQuadTetra::SMDS_DownQuadTetra(SMDS_UnstructuredGrid *grid) :
1217   SMDS_Down3D(grid, 4)
1218 {
1219   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1220   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1221   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1222   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1223 }
1224
1225 SMDS_DownQuadTetra::~SMDS_DownQuadTetra()
1226 {
1227 }
1228
1229 void SMDS_DownQuadTetra::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1230 {
1231   set<int> setNodes;
1232   setNodes.clear();
1233   for (int i = 0; i < orderedNodes.size(); i++)
1234     setNodes.insert(orderedNodes[i]);
1235   //MESSAGE("cellId = " << cellId);
1236
1237   vtkIdType npts = 0;
1238   vtkIdType *nodes; // will refer to the point id's of the volume
1239   _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
1240
1241   set<int> tofind;
1242   int ids[24] = { 0, 1, 2, 4, 5, 6,  0, 3, 1, 7, 8, 4,  2, 3, 0, 9, 7, 6,  1, 3, 2, 8, 9, 5 };
1243 //int ids[24] = { 2, 1, 0, 5, 4, 6,  1, 3, 0, 8, 7, 4,  0, 3, 2, 7, 9, 6,  2, 3, 1, 9, 8, 5 };
1244   for (int k = 0; k < 4; k++)
1245     {
1246       tofind.clear();
1247       for (int i = 0; i < 6; i++)
1248         tofind.insert(nodes[ids[6 * k + i]]);
1249       if (setNodes == tofind)
1250         {
1251           for (int i = 0; i < 6; i++)
1252             orderedNodes[i] = nodes[ids[6 * k + i]];
1253           return;
1254         }
1255     }
1256   MESSAGE("=== Problem volume " << _vtkCellIds[cellId] << " " << _grid->_mesh->fromVtkToSmds(_vtkCellIds[cellId]));
1257   MESSAGE(orderedNodes[0] << " " << orderedNodes[1] << " " << orderedNodes[2]);
1258   MESSAGE(nodes[0] << " " << nodes[1] << " " << nodes[2] << " " << nodes[3]);
1259 }
1260
1261 void SMDS_DownQuadTetra::addDownCell(int cellId, int lowCellId, unsigned char aType)
1262 {
1263   //ASSERT((cellId >=0)&& (cellId < _maxId));
1264   //ASSERT(aType == VTK_QUADRATIC_TRIANGLE);
1265   int *faces = &_cellIds[_nbDownCells * cellId];
1266   for (int i = 0; i < _nbDownCells; i++)
1267     {
1268       if (faces[i] < 0)
1269         {
1270           faces[i] = lowCellId;
1271           return;
1272         }
1273       if (faces[i] == lowCellId)
1274         return;
1275     }
1276   ASSERT(0);
1277 }
1278
1279 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1280  * The ordering of the ten points defining the quadratic tetrahedron cell is point id's (0-3,4-9)
1281  * where id's 0-3 are the four tetrahedron vertices;
1282  * and point id's 4-9 are the mid-edge nodes between (0,1), (1,2), (2,0), (0,3), (1,3), and (2,3).
1283  * @see vtkQuadraticTetra.h in Filtering.
1284  * @param cellId volumeId in vtkUnstructuredGrid
1285  * @param facesWithNodes vector of face descriptors to be filled
1286  */
1287 void SMDS_DownQuadTetra::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1288 {
1289   // --- find point id's of the volume
1290
1291   vtkIdType npts = 0;
1292   vtkIdType *nodes; // will refer to the point id's of the volume
1293   _grid->GetCellPoints(cellId, npts, nodes);
1294
1295   // --- create all the ordered list of node id's for each face
1296
1297   facesWithNodes.nbElems = 4;
1298
1299   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1300   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
1301   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
1302   facesWithNodes.elems[0].nodeIds[3] = nodes[4];
1303   facesWithNodes.elems[0].nodeIds[4] = nodes[5];
1304   facesWithNodes.elems[0].nodeIds[5] = nodes[6];
1305   facesWithNodes.elems[0].nbNodes = 6;
1306   facesWithNodes.elems[0].vtkType = VTK_QUADRATIC_TRIANGLE;
1307
1308   facesWithNodes.elems[1].nodeIds[0] = nodes[0];
1309   facesWithNodes.elems[1].nodeIds[1] = nodes[1];
1310   facesWithNodes.elems[1].nodeIds[2] = nodes[3];
1311   facesWithNodes.elems[1].nodeIds[3] = nodes[4];
1312   facesWithNodes.elems[1].nodeIds[4] = nodes[8];
1313   facesWithNodes.elems[1].nodeIds[5] = nodes[7];
1314   facesWithNodes.elems[1].nbNodes = 6;
1315   facesWithNodes.elems[1].vtkType = VTK_QUADRATIC_TRIANGLE;
1316
1317   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
1318   facesWithNodes.elems[2].nodeIds[1] = nodes[2];
1319   facesWithNodes.elems[2].nodeIds[2] = nodes[3];
1320   facesWithNodes.elems[2].nodeIds[3] = nodes[6];
1321   facesWithNodes.elems[2].nodeIds[4] = nodes[9];
1322   facesWithNodes.elems[2].nodeIds[5] = nodes[7];
1323   facesWithNodes.elems[2].nbNodes = 6;
1324   facesWithNodes.elems[2].vtkType = VTK_QUADRATIC_TRIANGLE;
1325
1326   facesWithNodes.elems[3].nodeIds[0] = nodes[1];
1327   facesWithNodes.elems[3].nodeIds[1] = nodes[2];
1328   facesWithNodes.elems[3].nodeIds[2] = nodes[3];
1329   facesWithNodes.elems[3].nodeIds[3] = nodes[5];
1330   facesWithNodes.elems[3].nodeIds[4] = nodes[9];
1331   facesWithNodes.elems[3].nodeIds[5] = nodes[8];
1332   facesWithNodes.elems[3].nbNodes = 6;
1333   facesWithNodes.elems[3].vtkType = VTK_QUADRATIC_TRIANGLE;
1334 }
1335
1336 // ---------------------------------------------------------------------------
1337
1338 SMDS_DownPyramid::SMDS_DownPyramid(SMDS_UnstructuredGrid *grid) :
1339   SMDS_Down3D(grid, 5)
1340 {
1341   _cellTypes.push_back(VTK_QUAD);
1342   _cellTypes.push_back(VTK_TRIANGLE);
1343   _cellTypes.push_back(VTK_TRIANGLE);
1344   _cellTypes.push_back(VTK_TRIANGLE);
1345   _cellTypes.push_back(VTK_TRIANGLE);
1346 }
1347
1348 SMDS_DownPyramid::~SMDS_DownPyramid()
1349 {
1350 }
1351
1352 void SMDS_DownPyramid::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1353 {
1354   // TODO
1355 }
1356
1357 void SMDS_DownPyramid::addDownCell(int cellId, int lowCellId, unsigned char aType)
1358 {
1359   //ASSERT((cellId >=0) && (cellId < _maxId));
1360   int *faces = &_cellIds[_nbDownCells * cellId];
1361   if (aType == VTK_QUAD)
1362     {
1363       if (faces[0] < 0)
1364         {
1365           faces[0] = lowCellId;
1366           return;
1367         }
1368       if (faces[0] == lowCellId)
1369         return;
1370     }
1371   else
1372     {
1373       //ASSERT(aType == VTK_TRIANGLE);
1374       for (int i = 1; i < _nbDownCells; i++)
1375         {
1376           if (faces[i] < 0)
1377             {
1378               faces[i] = lowCellId;
1379               return;
1380             }
1381           if (faces[i] == lowCellId)
1382             return;
1383         }
1384     }
1385   ASSERT(0);
1386 }
1387
1388 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1389  * The pyramid is defined by the five points (0-4) where (0,1,2,3) is the base of the pyramid which,
1390  * using the right hand rule, forms a quadrilateral whose normal points in the direction of the
1391  * pyramid apex at vertex #4.
1392  * @see vtkPyramid.h in Filtering.
1393  * @param cellId volumeId in vtkUnstructuredGrid
1394  * @param facesWithNodes vector of face descriptors to be filled
1395  */
1396 void SMDS_DownPyramid::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1397 {
1398   // --- find point id's of the volume
1399
1400   vtkIdType npts = 0;
1401   vtkIdType *nodes; // will refer to the point id's of the volume
1402   _grid->GetCellPoints(cellId, npts, nodes);
1403
1404   // --- create all the ordered list of node id's for each face
1405
1406   facesWithNodes.nbElems = 5;
1407
1408   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1409   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
1410   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
1411   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
1412   facesWithNodes.elems[0].nbNodes = 4;
1413   facesWithNodes.elems[0].vtkType = VTK_QUAD;
1414
1415   facesWithNodes.elems[1].nodeIds[0] = nodes[0];
1416   facesWithNodes.elems[1].nodeIds[1] = nodes[1];
1417   facesWithNodes.elems[1].nodeIds[2] = nodes[4];
1418   facesWithNodes.elems[1].nbNodes = 3;
1419   facesWithNodes.elems[1].vtkType = VTK_TRIANGLE;
1420
1421   facesWithNodes.elems[2].nodeIds[0] = nodes[1];
1422   facesWithNodes.elems[2].nodeIds[1] = nodes[2];
1423   facesWithNodes.elems[2].nodeIds[2] = nodes[4];
1424   facesWithNodes.elems[2].nbNodes = 3;
1425   facesWithNodes.elems[2].vtkType = VTK_TRIANGLE;
1426
1427   facesWithNodes.elems[3].nodeIds[0] = nodes[2];
1428   facesWithNodes.elems[3].nodeIds[1] = nodes[3];
1429   facesWithNodes.elems[3].nodeIds[2] = nodes[4];
1430   facesWithNodes.elems[3].nbNodes = 3;
1431   facesWithNodes.elems[3].vtkType = VTK_TRIANGLE;
1432
1433   facesWithNodes.elems[4].nodeIds[0] = nodes[3];
1434   facesWithNodes.elems[4].nodeIds[1] = nodes[0];
1435   facesWithNodes.elems[4].nodeIds[2] = nodes[4];
1436   facesWithNodes.elems[4].nbNodes = 3;
1437   facesWithNodes.elems[4].vtkType = VTK_TRIANGLE;
1438 }
1439
1440 // ---------------------------------------------------------------------------
1441
1442 SMDS_DownQuadPyramid::SMDS_DownQuadPyramid(SMDS_UnstructuredGrid *grid) :
1443   SMDS_Down3D(grid, 5)
1444 {
1445   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1446   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1447   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1448   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1449   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1450 }
1451
1452 SMDS_DownQuadPyramid::~SMDS_DownQuadPyramid()
1453 {
1454 }
1455
1456 void SMDS_DownQuadPyramid::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1457 {
1458   // TODO
1459 }
1460
1461 void SMDS_DownQuadPyramid::addDownCell(int cellId, int lowCellId, unsigned char aType)
1462 {
1463   //ASSERT((cellId >=0) && (cellId < _maxId));
1464   int *faces = &_cellIds[_nbDownCells * cellId];
1465   if (aType == VTK_QUADRATIC_QUAD)
1466     {
1467       if (faces[0] < 0)
1468         {
1469           faces[0] = lowCellId;
1470           return;
1471         }
1472       if (faces[0] == lowCellId)
1473         return;
1474     }
1475   else
1476     {
1477       //ASSERT(aType == VTK_QUADRATIC_TRIANGLE);
1478       for (int i = 1; i < _nbDownCells; i++)
1479         {
1480           if (faces[i] < 0)
1481             {
1482               faces[i] = lowCellId;
1483               return;
1484             }
1485           if (faces[i] == lowCellId)
1486             return;
1487         }
1488     }
1489   ASSERT(0);
1490 }
1491
1492 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1493  * The ordering of the thirteen points defining the quadratic pyramid cell is point id's (0-4,5-12)
1494  * where point id's 0-4 are the five corner vertices of the pyramid; followed
1495  * by eight mid-edge nodes (5-12). Note that these mid-edge nodes lie on the edges defined by
1496  * 5(0,1), 6(1,2), 7(2,3), 8(3,0), 9(0,4), 10(1,4), 11(2,4), 12(3,4).
1497  * @see vtkQuadraticPyramid.h in Filtering.
1498  * @param cellId volumeId in vtkUnstructuredGrid
1499  * @param facesWithNodes vector of face descriptors to be filled
1500  */
1501 void SMDS_DownQuadPyramid::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1502 {
1503   // --- find point id's of the volume
1504
1505   vtkIdType npts = 0;
1506   vtkIdType *nodes; // will refer to the point id's of the volume
1507   _grid->GetCellPoints(cellId, npts, nodes);
1508
1509   // --- create all the ordered list of node id's for each face
1510
1511   facesWithNodes.nbElems = 5;
1512
1513   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1514   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
1515   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
1516   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
1517   facesWithNodes.elems[0].nodeIds[4] = nodes[5];
1518   facesWithNodes.elems[0].nodeIds[5] = nodes[6];
1519   facesWithNodes.elems[0].nodeIds[6] = nodes[7];
1520   facesWithNodes.elems[0].nodeIds[7] = nodes[8];
1521   facesWithNodes.elems[0].nbNodes = 8;
1522   facesWithNodes.elems[0].vtkType = VTK_QUADRATIC_QUAD;
1523
1524   facesWithNodes.elems[1].nodeIds[0] = nodes[0];
1525   facesWithNodes.elems[1].nodeIds[1] = nodes[1];
1526   facesWithNodes.elems[1].nodeIds[2] = nodes[4];
1527   facesWithNodes.elems[1].nodeIds[3] = nodes[5];
1528   facesWithNodes.elems[1].nodeIds[4] = nodes[9];
1529   facesWithNodes.elems[1].nodeIds[5] = nodes[10];
1530   facesWithNodes.elems[1].nbNodes = 6;
1531   facesWithNodes.elems[1].vtkType = VTK_QUADRATIC_TRIANGLE;
1532
1533   facesWithNodes.elems[2].nodeIds[0] = nodes[1];
1534   facesWithNodes.elems[2].nodeIds[1] = nodes[2];
1535   facesWithNodes.elems[2].nodeIds[2] = nodes[4];
1536   facesWithNodes.elems[2].nodeIds[3] = nodes[6];
1537   facesWithNodes.elems[2].nodeIds[4] = nodes[10];
1538   facesWithNodes.elems[2].nodeIds[5] = nodes[11];
1539   facesWithNodes.elems[2].nbNodes = 6;
1540   facesWithNodes.elems[2].vtkType = VTK_QUADRATIC_TRIANGLE;
1541
1542   facesWithNodes.elems[3].nodeIds[0] = nodes[2];
1543   facesWithNodes.elems[3].nodeIds[1] = nodes[3];
1544   facesWithNodes.elems[3].nodeIds[2] = nodes[4];
1545   facesWithNodes.elems[3].nodeIds[3] = nodes[7];
1546   facesWithNodes.elems[3].nodeIds[4] = nodes[11];
1547   facesWithNodes.elems[3].nodeIds[5] = nodes[12];
1548   facesWithNodes.elems[3].nbNodes = 6;
1549   facesWithNodes.elems[3].vtkType = VTK_QUADRATIC_TRIANGLE;
1550
1551   facesWithNodes.elems[4].nodeIds[0] = nodes[3];
1552   facesWithNodes.elems[4].nodeIds[1] = nodes[0];
1553   facesWithNodes.elems[4].nodeIds[2] = nodes[4];
1554   facesWithNodes.elems[4].nodeIds[3] = nodes[8];
1555   facesWithNodes.elems[4].nodeIds[4] = nodes[9];
1556   facesWithNodes.elems[4].nodeIds[5] = nodes[12];
1557   facesWithNodes.elems[4].nbNodes = 6;
1558   facesWithNodes.elems[4].vtkType = VTK_QUADRATIC_TRIANGLE;
1559 }
1560
1561 // ---------------------------------------------------------------------------
1562
1563 SMDS_DownPenta::SMDS_DownPenta(SMDS_UnstructuredGrid *grid) :
1564   SMDS_Down3D(grid, 5)
1565 {
1566   _cellTypes.push_back(VTK_QUAD);
1567   _cellTypes.push_back(VTK_QUAD);
1568   _cellTypes.push_back(VTK_QUAD);
1569   _cellTypes.push_back(VTK_TRIANGLE);
1570   _cellTypes.push_back(VTK_TRIANGLE);
1571 }
1572
1573 SMDS_DownPenta::~SMDS_DownPenta()
1574 {
1575 }
1576
1577 void SMDS_DownPenta::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1578 {
1579   // TODO
1580 }
1581
1582 void SMDS_DownPenta::addDownCell(int cellId, int lowCellId, unsigned char aType)
1583 {
1584   //ASSERT((cellId >=0) && (cellId < _maxId));
1585   int *faces = &_cellIds[_nbDownCells * cellId];
1586   if (aType == VTK_QUAD)
1587     for (int i = 0; i < 2; i++)
1588       {
1589         if (faces[i] < 0)
1590           {
1591             faces[i] = lowCellId;
1592             return;
1593           }
1594         if (faces[i] == lowCellId)
1595           return;
1596       }
1597   else
1598     {
1599       //ASSERT(aType == VTK_TRIANGLE);
1600       for (int i = 2; i < _nbDownCells; i++)
1601         {
1602           if (faces[i] < 0)
1603             {
1604               faces[i] = lowCellId;
1605               return;
1606             }
1607           if (faces[i] == lowCellId)
1608             return;
1609         }
1610     }
1611   ASSERT(0);
1612 }
1613
1614 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's.
1615  * A wedge or pentahedron consists of two triangular and three quadrilateral faces
1616  * and is defined by the six points (0-5) where (0,1,2) is the base of the wedge which,
1617  * using the right hand rule, forms a triangle whose normal points outward
1618  * (away from the triangular face (3,4,5)).
1619  * @see vtkWedge.h in Filtering
1620  * @param cellId volumeId in vtkUnstructuredGrid
1621  * @param facesWithNodes vector of face descriptors to be filled
1622  */
1623 void SMDS_DownPenta::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1624 {
1625   // --- find point id's of the volume
1626
1627   vtkIdType npts = 0;
1628   vtkIdType *nodes; // will refer to the point id's of the volume
1629   _grid->GetCellPoints(cellId, npts, nodes);
1630
1631   // --- create all the ordered list of node id's for each face
1632
1633   facesWithNodes.nbElems = 5;
1634
1635   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1636   facesWithNodes.elems[0].nodeIds[1] = nodes[2];
1637   facesWithNodes.elems[0].nodeIds[2] = nodes[5];
1638   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
1639   facesWithNodes.elems[0].nbNodes = 4;
1640   facesWithNodes.elems[0].vtkType = VTK_QUAD;
1641
1642   facesWithNodes.elems[1].nodeIds[0] = nodes[1];
1643   facesWithNodes.elems[1].nodeIds[1] = nodes[2];
1644   facesWithNodes.elems[1].nodeIds[2] = nodes[5];
1645   facesWithNodes.elems[1].nodeIds[3] = nodes[4];
1646   facesWithNodes.elems[1].nbNodes = 4;
1647   facesWithNodes.elems[1].vtkType = VTK_QUAD;
1648
1649   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
1650   facesWithNodes.elems[2].nodeIds[1] = nodes[1];
1651   facesWithNodes.elems[2].nodeIds[2] = nodes[4];
1652   facesWithNodes.elems[2].nodeIds[3] = nodes[3];
1653   facesWithNodes.elems[2].nbNodes = 4;
1654   facesWithNodes.elems[2].vtkType = VTK_QUAD;
1655
1656   facesWithNodes.elems[3].nodeIds[0] = nodes[0];
1657   facesWithNodes.elems[3].nodeIds[1] = nodes[1];
1658   facesWithNodes.elems[3].nodeIds[2] = nodes[2];
1659   facesWithNodes.elems[3].nbNodes = 3;
1660   facesWithNodes.elems[3].vtkType = VTK_TRIANGLE;
1661
1662   facesWithNodes.elems[4].nodeIds[0] = nodes[3];
1663   facesWithNodes.elems[4].nodeIds[1] = nodes[4];
1664   facesWithNodes.elems[4].nodeIds[2] = nodes[5];
1665   facesWithNodes.elems[4].nbNodes = 3;
1666   facesWithNodes.elems[4].vtkType = VTK_TRIANGLE;
1667 }
1668
1669 // ---------------------------------------------------------------------------
1670
1671 SMDS_DownQuadPenta::SMDS_DownQuadPenta(SMDS_UnstructuredGrid *grid) :
1672   SMDS_Down3D(grid, 5)
1673 {
1674   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1675   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1676   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1677   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1678   _cellTypes.push_back(VTK_QUADRATIC_TRIANGLE);
1679 }
1680
1681 SMDS_DownQuadPenta::~SMDS_DownQuadPenta()
1682 {
1683 }
1684
1685 void SMDS_DownQuadPenta::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1686 {
1687   // TODO
1688 }
1689
1690 void SMDS_DownQuadPenta::addDownCell(int cellId, int lowCellId, unsigned char aType)
1691 {
1692   //ASSERT((cellId >=0) && (cellId < _maxId));
1693   int *faces = &_cellIds[_nbDownCells * cellId];
1694   if (aType == VTK_QUADRATIC_QUAD)
1695     for (int i = 0; i < 2; i++)
1696       {
1697         if (faces[i] < 0)
1698           {
1699             faces[i] = lowCellId;
1700             return;
1701           }
1702         if (faces[i] == lowCellId)
1703           return;
1704       }
1705   else
1706     {
1707       //ASSERT(aType == VTK_QUADRATIC_TRIANGLE);
1708       for (int i = 2; i < _nbDownCells; i++)
1709         {
1710           if (faces[i] < 0)
1711             {
1712               faces[i] = lowCellId;
1713               return;
1714             }
1715           if (faces[i] == lowCellId)
1716             return;
1717         }
1718     }
1719   ASSERT(0);
1720 }
1721
1722 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1723  * The quadratic wedge (or pentahedron) is defined by fifteen points.
1724  * The ordering of the fifteen points defining the cell is point id's (0-5,6-14)
1725  * where point id's 0-5 are the six corner vertices of the wedge, followed by
1726  * nine mid-edge nodes (6-14). Note that these mid-edge nodes lie on the edges defined by
1727  * (0,1), (1,2), (2,0), (3,4), (4,5), (5,3), (0,3), (1,4), (2,5).
1728  * @see vtkQuadraticWedge.h in Filtering
1729  * @param cellId volumeId in vtkUnstructuredGrid
1730  * @param facesWithNodes vector of face descriptors to be filled
1731  */
1732 void SMDS_DownQuadPenta::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1733 {
1734   // --- find point id's of the volume
1735
1736   vtkIdType npts = 0;
1737   vtkIdType *nodes; // will refer to the point id's of the volume
1738   _grid->GetCellPoints(cellId, npts, nodes);
1739
1740   // --- create all the ordered list of node id's for each face
1741
1742   facesWithNodes.nbElems = 5;
1743
1744   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1745   facesWithNodes.elems[0].nodeIds[1] = nodes[2];
1746   facesWithNodes.elems[0].nodeIds[2] = nodes[5];
1747   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
1748   facesWithNodes.elems[0].nodeIds[4] = nodes[8];
1749   facesWithNodes.elems[0].nodeIds[5] = nodes[14];
1750   facesWithNodes.elems[0].nodeIds[6] = nodes[11];
1751   facesWithNodes.elems[0].nodeIds[7] = nodes[12];
1752   facesWithNodes.elems[0].nbNodes = 8;
1753   facesWithNodes.elems[0].vtkType = VTK_QUADRATIC_QUAD;
1754
1755   facesWithNodes.elems[1].nodeIds[0] = nodes[1];
1756   facesWithNodes.elems[1].nodeIds[1] = nodes[2];
1757   facesWithNodes.elems[1].nodeIds[2] = nodes[5];
1758   facesWithNodes.elems[1].nodeIds[3] = nodes[4];
1759   facesWithNodes.elems[1].nodeIds[4] = nodes[7];
1760   facesWithNodes.elems[1].nodeIds[5] = nodes[14];
1761   facesWithNodes.elems[1].nodeIds[6] = nodes[10];
1762   facesWithNodes.elems[1].nodeIds[7] = nodes[13];
1763   facesWithNodes.elems[1].nbNodes = 8;
1764   facesWithNodes.elems[1].vtkType = VTK_QUADRATIC_QUAD;
1765
1766   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
1767   facesWithNodes.elems[2].nodeIds[1] = nodes[1];
1768   facesWithNodes.elems[2].nodeIds[2] = nodes[4];
1769   facesWithNodes.elems[2].nodeIds[3] = nodes[3];
1770   facesWithNodes.elems[2].nodeIds[4] = nodes[6];
1771   facesWithNodes.elems[2].nodeIds[5] = nodes[13];
1772   facesWithNodes.elems[2].nodeIds[6] = nodes[9];
1773   facesWithNodes.elems[2].nodeIds[7] = nodes[12];
1774   facesWithNodes.elems[2].nbNodes = 8;
1775   facesWithNodes.elems[2].vtkType = VTK_QUADRATIC_QUAD;
1776
1777   facesWithNodes.elems[3].nodeIds[0] = nodes[0];
1778   facesWithNodes.elems[3].nodeIds[1] = nodes[1];
1779   facesWithNodes.elems[3].nodeIds[2] = nodes[2];
1780   facesWithNodes.elems[3].nodeIds[3] = nodes[6];
1781   facesWithNodes.elems[3].nodeIds[4] = nodes[7];
1782   facesWithNodes.elems[3].nodeIds[5] = nodes[8];
1783   facesWithNodes.elems[3].nbNodes = 6;
1784   facesWithNodes.elems[3].vtkType = VTK_QUADRATIC_TRIANGLE;
1785
1786   facesWithNodes.elems[4].nodeIds[0] = nodes[3];
1787   facesWithNodes.elems[4].nodeIds[1] = nodes[4];
1788   facesWithNodes.elems[4].nodeIds[2] = nodes[5];
1789   facesWithNodes.elems[4].nodeIds[3] = nodes[9];
1790   facesWithNodes.elems[4].nodeIds[4] = nodes[10];
1791   facesWithNodes.elems[4].nodeIds[5] = nodes[11];
1792   facesWithNodes.elems[4].nbNodes = 6;
1793   facesWithNodes.elems[4].vtkType = VTK_QUADRATIC_TRIANGLE;
1794 }
1795
1796 // ---------------------------------------------------------------------------
1797
1798 SMDS_DownHexa::SMDS_DownHexa(SMDS_UnstructuredGrid *grid) :
1799   SMDS_Down3D(grid, 6)
1800 {
1801   _cellTypes.push_back(VTK_QUAD);
1802   _cellTypes.push_back(VTK_QUAD);
1803   _cellTypes.push_back(VTK_QUAD);
1804   _cellTypes.push_back(VTK_QUAD);
1805   _cellTypes.push_back(VTK_QUAD);
1806   _cellTypes.push_back(VTK_QUAD);
1807 }
1808
1809 SMDS_DownHexa::~SMDS_DownHexa()
1810 {
1811 }
1812
1813 void SMDS_DownHexa::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1814 {
1815   set<int> setNodes;
1816   setNodes.clear();
1817   for (int i = 0; i < orderedNodes.size(); i++)
1818     setNodes.insert(orderedNodes[i]);
1819   //MESSAGE("cellId = " << cellId);
1820
1821   vtkIdType npts = 0;
1822   vtkIdType *nodes; // will refer to the point id's of the volume
1823   _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
1824
1825   set<int> tofind;
1826 //int ids[24] = { 0, 1, 2, 3,  7, 6, 5, 4,  4, 0, 3, 7,  5, 1, 0, 4,  6, 2, 1, 5,  7, 3, 2, 6};
1827   int ids[24] = { 3, 2, 1, 0,  4, 5, 6, 7,  7, 3, 0, 4,  4, 0, 1, 5,  5, 1, 2, 6,  6, 2, 3, 7};
1828   for (int k = 0; k < 6; k++) // loop on the 6 faces
1829     {
1830       tofind.clear();
1831       for (int i = 0; i < 4; i++)
1832         tofind.insert(nodes[ids[4 * k + i]]); // node ids of the face i
1833       if (setNodes == tofind)
1834         {
1835           for (int i = 0; i < 4; i++)
1836             orderedNodes[i] = nodes[ids[4 * k + i]];
1837           return;
1838         }
1839     }
1840   MESSAGE("=== Problem volume " << _vtkCellIds[cellId] << " " << _grid->_mesh->fromVtkToSmds(_vtkCellIds[cellId]));
1841   MESSAGE(orderedNodes[0] << " " << orderedNodes[1] << " " << orderedNodes[2] << " " << orderedNodes[3]);
1842   MESSAGE(nodes[0] << " " << nodes[1] << " " << nodes[2] << " " << nodes[3]);
1843   MESSAGE(nodes[4] << " " << nodes[5] << " " << nodes[6] << " " << nodes[7]);
1844 }
1845
1846 void SMDS_DownHexa::addDownCell(int cellId, int lowCellId, unsigned char aType)
1847 {
1848   //ASSERT((cellId >=0)&& (cellId < _maxId));
1849   int *faces = &_cellIds[_nbDownCells * cellId];
1850   for (int i = 0; i < _nbDownCells; i++)
1851     {
1852       if (faces[i] < 0)
1853         {
1854           faces[i] = lowCellId;
1855           return;
1856         }
1857       if (faces[i] == lowCellId)
1858         return;
1859     }
1860   ASSERT(0);
1861   // MESSAGE("-------------------------------------> trop de faces ! " << cellId << " " << lowCellId);
1862 }
1863
1864 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1865  * The hexahedron is defined by the eight points (0-7), where (0,1,2,3) is the base
1866  * of the hexahedron which, using the right hand rule, forms a quadrilateral whose normal
1867  * points in the direction of the opposite face (4,5,6,7).
1868  * @see vtkHexahedron.h in Filtering
1869  * @param cellId volumeId in vtkUnstructuredGrid
1870  * @param facesWithNodes vector of face descriptors to be filled
1871  */
1872 void SMDS_DownHexa::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1873 {
1874   // --- find point id's of the volume
1875
1876   vtkIdType npts = 0;
1877   vtkIdType *nodes; // will refer to the point id's of the volume
1878   _grid->GetCellPoints(cellId, npts, nodes);
1879
1880   // --- create all the ordered list of node id's for each face
1881
1882   facesWithNodes.nbElems = 6;
1883
1884   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
1885   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
1886   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
1887   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
1888   facesWithNodes.elems[0].nbNodes = 4;
1889   facesWithNodes.elems[0].vtkType = VTK_QUAD;
1890
1891   facesWithNodes.elems[1].nodeIds[0] = nodes[4];
1892   facesWithNodes.elems[1].nodeIds[1] = nodes[5];
1893   facesWithNodes.elems[1].nodeIds[2] = nodes[6];
1894   facesWithNodes.elems[1].nodeIds[3] = nodes[7];
1895   facesWithNodes.elems[1].nbNodes = 4;
1896   facesWithNodes.elems[1].vtkType = VTK_QUAD;
1897
1898   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
1899   facesWithNodes.elems[2].nodeIds[1] = nodes[1];
1900   facesWithNodes.elems[2].nodeIds[2] = nodes[5];
1901   facesWithNodes.elems[2].nodeIds[3] = nodes[4];
1902   facesWithNodes.elems[2].nbNodes = 4;
1903   facesWithNodes.elems[2].vtkType = VTK_QUAD;
1904
1905   facesWithNodes.elems[3].nodeIds[0] = nodes[1];
1906   facesWithNodes.elems[3].nodeIds[1] = nodes[2];
1907   facesWithNodes.elems[3].nodeIds[2] = nodes[6];
1908   facesWithNodes.elems[3].nodeIds[3] = nodes[5];
1909   facesWithNodes.elems[3].nbNodes = 4;
1910   facesWithNodes.elems[3].vtkType = VTK_QUAD;
1911
1912   facesWithNodes.elems[4].nodeIds[0] = nodes[2];
1913   facesWithNodes.elems[4].nodeIds[1] = nodes[6];
1914   facesWithNodes.elems[4].nodeIds[2] = nodes[7];
1915   facesWithNodes.elems[4].nodeIds[3] = nodes[3];
1916   facesWithNodes.elems[4].nbNodes = 4;
1917   facesWithNodes.elems[4].vtkType = VTK_QUAD;
1918
1919   facesWithNodes.elems[5].nodeIds[0] = nodes[3];
1920   facesWithNodes.elems[5].nodeIds[1] = nodes[7];
1921   facesWithNodes.elems[5].nodeIds[2] = nodes[4];
1922   facesWithNodes.elems[5].nodeIds[3] = nodes[0];
1923   facesWithNodes.elems[5].nbNodes = 4;
1924   facesWithNodes.elems[5].vtkType = VTK_QUAD;
1925 }
1926
1927 // ---------------------------------------------------------------------------
1928
1929 SMDS_DownQuadHexa::SMDS_DownQuadHexa(SMDS_UnstructuredGrid *grid) :
1930   SMDS_Down3D(grid, 6)
1931 {
1932   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1933   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1934   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1935   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1936   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1937   _cellTypes.push_back(VTK_QUADRATIC_QUAD);
1938 }
1939
1940 SMDS_DownQuadHexa::~SMDS_DownQuadHexa()
1941 {
1942 }
1943
1944 void SMDS_DownQuadHexa::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1945 {
1946   set<int> setNodes;
1947   setNodes.clear();
1948   for (int i = 0; i < orderedNodes.size(); i++)
1949     setNodes.insert(orderedNodes[i]);
1950   //MESSAGE("cellId = " << cellId);
1951
1952   vtkIdType npts = 0;
1953   vtkIdType *nodes; // will refer to the point id's of the volume
1954   _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
1955
1956   set<int> tofind;
1957   //int ids[24] = { 3, 2, 1, 0,  4, 5, 6, 7,  7, 3, 0, 4,  4, 0, 1, 5,  5, 1, 2, 6,  6, 2, 3, 7};
1958   int ids[48] = { 3, 2, 1, 0,10, 9, 8,11,   4, 5, 6, 7,12,13,14,15,   7, 3, 0, 4,19,11,16,15,
1959                   4, 0, 1, 5,16, 8,17,12,   5, 1, 2, 6,17, 9,18,13,   6, 2, 3, 7,18,10,19,14};
1960   for (int k = 0; k < 6; k++)
1961     {
1962       tofind.clear();
1963       for (int i = 0; i < 8; i++)
1964         tofind.insert(nodes[ids[8 * k + i]]);
1965       if (setNodes == tofind)
1966         {
1967           for (int i = 0; i < 8; i++)
1968             orderedNodes[i] = nodes[ids[8 * k + i]];
1969           return;
1970         }
1971     }
1972   MESSAGE("=== Problem volume " << _vtkCellIds[cellId] << " " << _grid->_mesh->fromVtkToSmds(_vtkCellIds[cellId]));
1973   MESSAGE(orderedNodes[0] << " " << orderedNodes[1] << " " << orderedNodes[2] << " " << orderedNodes[3]);
1974   MESSAGE(nodes[0] << " " << nodes[1] << " " << nodes[2] << " " << nodes[3]);
1975 }
1976
1977 void SMDS_DownQuadHexa::addDownCell(int cellId, int lowCellId, unsigned char aType)
1978 {
1979   //ASSERT((cellId >=0)&& (cellId < _maxId));
1980   int *faces = &_cellIds[_nbDownCells * cellId];
1981   for (int i = 0; i < _nbDownCells; i++)
1982     {
1983       if (faces[i] < 0)
1984         {
1985           faces[i] = lowCellId;
1986           return;
1987         }
1988       if (faces[i] == lowCellId)
1989         return;
1990     }
1991   ASSERT(0);
1992 }
1993
1994 /*! Create a list of faces described by a vtk Type and  an ordered set of Node Id's
1995  * The ordering of the twenty points defining the quadratic hexahedron cell is point id's (0-7,8-19)
1996  * where point id's 0-7 are the eight corner vertices of the cube, followed by twelve mid-edge nodes (8-19).
1997  * Note that these mid-edge nodes lie on the edges defined by
1998  * (0,1), (1,2), (2,3), (3,0), (4,5), (5,6), (6,7), (7,4), (0,4), (1,5), (2,6), (3,7).
1999  * @see vtkQuadraticHexahedron.h in Filtering
2000  * @param cellId volumeId in vtkUnstructuredGrid
2001  * @param facesWithNodes vector of face descriptors to be filled
2002  */
2003 void SMDS_DownQuadHexa::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
2004 {
2005   // --- find point id's of the volume
2006
2007   vtkIdType npts = 0;
2008   vtkIdType *nodes; // will refer to the point id's of the volume
2009   _grid->GetCellPoints(cellId, npts, nodes);
2010
2011   // --- create all the ordered list of node id's for each face
2012
2013   facesWithNodes.nbElems = 6;
2014
2015   facesWithNodes.elems[0].nodeIds[0] = nodes[0];
2016   facesWithNodes.elems[0].nodeIds[1] = nodes[1];
2017   facesWithNodes.elems[0].nodeIds[2] = nodes[2];
2018   facesWithNodes.elems[0].nodeIds[3] = nodes[3];
2019   facesWithNodes.elems[0].nodeIds[4] = nodes[8];
2020   facesWithNodes.elems[0].nodeIds[5] = nodes[9];
2021   facesWithNodes.elems[0].nodeIds[6] = nodes[10];
2022   facesWithNodes.elems[0].nodeIds[7] = nodes[11];
2023   facesWithNodes.elems[0].nbNodes = 8;
2024   facesWithNodes.elems[0].vtkType = VTK_QUADRATIC_QUAD;
2025
2026   facesWithNodes.elems[1].nodeIds[0] = nodes[4];
2027   facesWithNodes.elems[1].nodeIds[1] = nodes[5];
2028   facesWithNodes.elems[1].nodeIds[2] = nodes[6];
2029   facesWithNodes.elems[1].nodeIds[3] = nodes[7];
2030   facesWithNodes.elems[1].nodeIds[4] = nodes[12];
2031   facesWithNodes.elems[1].nodeIds[5] = nodes[13];
2032   facesWithNodes.elems[1].nodeIds[6] = nodes[14];
2033   facesWithNodes.elems[1].nodeIds[7] = nodes[15];
2034   facesWithNodes.elems[1].nbNodes = 8;
2035   facesWithNodes.elems[1].vtkType = VTK_QUADRATIC_QUAD;
2036
2037   facesWithNodes.elems[2].nodeIds[0] = nodes[0];
2038   facesWithNodes.elems[2].nodeIds[1] = nodes[1];
2039   facesWithNodes.elems[2].nodeIds[2] = nodes[5];
2040   facesWithNodes.elems[2].nodeIds[3] = nodes[4];
2041   facesWithNodes.elems[2].nodeIds[4] = nodes[8];
2042   facesWithNodes.elems[2].nodeIds[5] = nodes[17];
2043   facesWithNodes.elems[2].nodeIds[6] = nodes[12];
2044   facesWithNodes.elems[2].nodeIds[7] = nodes[16];
2045   facesWithNodes.elems[2].nbNodes = 8;
2046   facesWithNodes.elems[2].vtkType = VTK_QUADRATIC_QUAD;
2047
2048   facesWithNodes.elems[3].nodeIds[0] = nodes[1];
2049   facesWithNodes.elems[3].nodeIds[1] = nodes[2];
2050   facesWithNodes.elems[3].nodeIds[2] = nodes[6];
2051   facesWithNodes.elems[3].nodeIds[3] = nodes[5];
2052   facesWithNodes.elems[3].nodeIds[4] = nodes[9];
2053   facesWithNodes.elems[3].nodeIds[5] = nodes[18];
2054   facesWithNodes.elems[3].nodeIds[6] = nodes[13];
2055   facesWithNodes.elems[3].nodeIds[7] = nodes[17];
2056   facesWithNodes.elems[3].nbNodes = 8;
2057   facesWithNodes.elems[3].vtkType = VTK_QUADRATIC_QUAD;
2058
2059   facesWithNodes.elems[4].nodeIds[0] = nodes[2];
2060   facesWithNodes.elems[4].nodeIds[1] = nodes[6];
2061   facesWithNodes.elems[4].nodeIds[2] = nodes[7];
2062   facesWithNodes.elems[4].nodeIds[3] = nodes[3];
2063   facesWithNodes.elems[4].nodeIds[4] = nodes[18];
2064   facesWithNodes.elems[4].nodeIds[5] = nodes[14];
2065   facesWithNodes.elems[4].nodeIds[6] = nodes[19];
2066   facesWithNodes.elems[4].nodeIds[7] = nodes[10];
2067   facesWithNodes.elems[4].nbNodes = 8;
2068   facesWithNodes.elems[4].vtkType = VTK_QUADRATIC_QUAD;
2069
2070   facesWithNodes.elems[5].nodeIds[0] = nodes[3];
2071   facesWithNodes.elems[5].nodeIds[1] = nodes[7];
2072   facesWithNodes.elems[5].nodeIds[2] = nodes[4];
2073   facesWithNodes.elems[5].nodeIds[3] = nodes[0];
2074   facesWithNodes.elems[5].nodeIds[4] = nodes[19];
2075   facesWithNodes.elems[5].nodeIds[5] = nodes[15];
2076   facesWithNodes.elems[5].nodeIds[6] = nodes[16];
2077   facesWithNodes.elems[5].nodeIds[7] = nodes[11];
2078   facesWithNodes.elems[5].nbNodes = 8;
2079   facesWithNodes.elems[5].vtkType = VTK_QUADRATIC_QUAD;
2080 }
2081
2082 // ---------------------------------------------------------------------------
2083