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