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