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