1 // Copyright (C) 2010-2016 CEA/DEN, EDF R&D, OPEN CASCADE
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 // File: SMDS_Downward.cxx
21 // Created: Jun 3, 2010
24 #include "SMDS_Downward.hxx"
25 #include "SMDS_Mesh.hxx"
26 #include "utilities.h"
28 #include <vtkCellType.h>
29 #include <vtkCellLinks.h>
35 // ---------------------------------------------------------------------------
37 vector<int> SMDS_Downward::_cellDimension;
39 /*! get the dimension of a cell (1,2,3 for 1D, 2D 3D) given the vtk cell type
41 * @param cellType vtk cell type @see vtkCellType.h
44 int SMDS_Downward::getCellDimension(unsigned char cellType)
46 if (_cellDimension.empty())
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;
69 return _cellDimension[cellType];
72 // ---------------------------------------------------------------------------
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.
80 SMDS_Downward::SMDS_Downward(SMDS_UnstructuredGrid *grid, int nbDownCells) :
81 _grid(grid), _nbDownCells(nbDownCells)
84 this->_cellIds.clear();
85 this->_cellTypes.clear();
86 if (_cellDimension.empty())
87 getCellDimension( VTK_LINE );
90 SMDS_Downward::~SMDS_Downward()
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.
101 int SMDS_Downward::addCell(int vtkId)
105 localId = _grid->CellIdToDownId(vtkId);
109 localId = this->_maxId;
111 this->allocate(_maxId);
114 this->_vtkCellIds[localId] = vtkId;
115 _grid->setCellIdToDownId(vtkId, localId);
117 this->initCell(localId);
121 /*! generic method do nothing. see derived methods
125 void SMDS_Downward::initCell(int cellId)
129 /*! Get the number of downward entities associated to a cell (always the same for a given vtk type of cell)
131 * @param cellId not used here.
134 int SMDS_Downward::getNumberOfDownCells(int cellId)
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.
145 const int* SMDS_Downward::getDownCells(int cellId)
147 //ASSERT((cellId >=0) && (cellId < _maxId));
148 return &_cellIds[_nbDownCells * cellId];
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).
154 * @param cellId index of the cell in the downward structure relative to a vtk cell type.
155 * @return table of downward entities types.
157 const unsigned char* SMDS_Downward::getDownTypes(int cellId)
159 return &_cellTypes[0];
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).
168 void SMDS_Downward::addDownCell(int cellId, int lowCellId, unsigned char aType)
170 ASSERT(0); // must be re-implemented in derived class
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).
179 void SMDS_Downward::addUpCell(int cellId, int upCellId, unsigned char aType)
181 ASSERT(0); // must be re-implemented in derived class
184 int SMDS_Downward::getNodeSet(int cellId, int* nodeSet)
189 // ---------------------------------------------------------------------------
191 SMDS_Down1D::SMDS_Down1D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
192 SMDS_Downward(grid, nbDownCells)
194 _upCellIdsVector.clear();
195 _upCellTypesVector.clear();
197 _upCellTypes.clear();
198 _upCellIndex.clear();
201 SMDS_Down1D::~SMDS_Down1D()
205 /*! clear vectors used to reference 2D cells containing the edge
209 void SMDS_Down1D::initCell(int cellId)
211 _upCellIdsVector[cellId].clear();
212 _upCellTypesVector[cellId].clear();
215 /*! Resize the downward connectivity storage vector if needed.
217 * @param nbElems total number of elements of the same type required
219 void SMDS_Down1D::allocate(int nbElems)
221 if (nbElems >= (int)_vtkCellIds.size())
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);
230 void SMDS_Down1D::compactStorage()
232 _cellIds.resize(_nbDownCells * _maxId);
233 _vtkCellIds.resize(_maxId);
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] [
242 for (int i = 0; i < _maxId; i++)
244 _upCellIndex[i] = current;
245 for (size_t j = 0; j < _upCellIdsVector[i].size(); j++)
247 _upCellIds[current] = _upCellIdsVector[i][j];
248 _upCellTypes[current] = _upCellTypesVector[i][j];
252 _upCellIndex[_maxId] = current;
254 _upCellIdsVector.clear();
255 _upCellTypesVector.clear();
258 void SMDS_Down1D::addUpCell(int cellId, int upCellId, unsigned char aType)
260 //ASSERT((cellId >=0) && (cellId < _maxId));
261 int nbFaces = _upCellIdsVector[cellId].size();
262 for (int i = 0; i < nbFaces; i++)
264 if ((_upCellIdsVector[cellId][i] == upCellId) && (_upCellTypesVector[cellId][i] == aType))
266 return; // already done
269 _upCellIdsVector[cellId].push_back(upCellId);
270 _upCellTypesVector[cellId].push_back(aType);
273 int SMDS_Down1D::getNumberOfUpCells(int cellId)
275 //ASSERT((cellId >=0) && (cellId < _maxId));
276 return _upCellIndex[cellId + 1] - _upCellIndex[cellId];
279 const int* SMDS_Down1D::getUpCells(int cellId)
281 //ASSERT((cellId >=0) && (cellId < _maxId));
282 return &_upCellIds[_upCellIndex[cellId]];
285 const unsigned char* SMDS_Down1D::getUpTypes(int cellId)
287 //ASSERT((cellId >=0) && (cellId < _maxId));
288 return &_upCellTypes[_upCellIndex[cellId]];
291 void SMDS_Down1D::getNodeIds(int cellId, std::set<int>& nodeSet)
293 for (int i = 0; i < _nbDownCells; i++)
294 nodeSet.insert(_cellIds[_nbDownCells * cellId + i]);
297 int SMDS_Down1D::getNodeSet(int cellId, int* nodeSet)
299 for (int i = 0; i < _nbDownCells; i++)
300 nodeSet[i] = _cellIds[_nbDownCells * cellId + i];
304 void SMDS_Down1D::setNodes(int cellId, int vtkId)
307 vtkIdType *pts; // 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++)
313 _cellIds[_nbDownCells * cellId + i] = pts[i];
317 void SMDS_Down1D::setNodes(int cellId, const int* nodeIds)
319 //ASSERT(nodeIds.size() == _nbDownCells);
320 for (int i = 0; i < _nbDownCells; i++)
322 _cellIds[_nbDownCells * cellId + i] = nodeIds[i];
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)
332 int SMDS_Down1D::computeVtkCells(int cellId, std::vector<int>& vtkIds)
336 // --- find all the cells the points belong to, and how many of the points belong to a given cell
338 int *pts = &_cellIds[_nbDownCells * cellId];
339 int ncells = this->computeVtkCells(pts, vtkIds);
343 /*! Build the list of vtkUnstructuredGrid cells containing the edge.
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)
349 int SMDS_Down1D::computeVtkCells(int *pts, std::vector<int>& vtkIds)
352 // --- find all the cells the points belong to, and how many of the points belong to a given cell
357 for (int i = 0; i < _nbDownCells; i++)
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++)
364 int vtkCellId = cells[j];
366 for (int k = 0; k < cnt; k++)
368 if (cellIds[k] == vtkCellId)
377 cellIds[cnt] = vtkCellId;
379 // TODO ASSERT(cnt<1000);
385 // --- find the face and volume cells: they contains all the points and are of type volume or face
388 for (int i = 0; i < cnt; i++)
390 if (cellCnt[i] == _nbDownCells)
392 int vtkElemId = cellIds[i];
393 int vtkType = _grid->GetCellType(vtkElemId);
394 if (SMDS_Downward::getCellDimension(vtkType) > 1)
396 vtkIds.push_back(vtkElemId);
405 /*! Build the list of downward faces from a list of vtk cells.
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
413 int SMDS_Down1D::computeFaces(int cellId, int* vtkIds, int nbcells, int* downFaces, unsigned char* downTypes)
415 int *pts = &_cellIds[_nbDownCells * cellId];
416 int nbFaces = this->computeFaces(pts, vtkIds, nbcells, downFaces, downTypes);
420 /*! Build the list of downward faces from a list of vtk cells.
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
428 int SMDS_Down1D::computeFaces(int* pts, int* vtkIds, int nbcells, int* downFaces, unsigned char* downTypes)
431 for (int i = 0; i < nbcells; i++)
433 int vtkId = vtkIds[i];
434 int vtkType = _grid->GetCellType(vtkId);
435 if (SMDS_Downward::getCellDimension(vtkType) == 2)
437 int faceId = _grid->CellIdToDownId(vtkId);
438 downFaces[cnt] = faceId;
439 downTypes[cnt] = vtkType;
442 else if (SMDS_Downward::getCellDimension(vtkType) == 3)
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++)
452 SMDS_Down2D *downFace = static_cast<SMDS_Down2D*> (_grid->getDownArray(downTypesVol[n]));
453 bool isInFace = downFace->isInFace(faceIds[n], pts, _nbDownCells);
456 bool alreadySet = false;
457 for (int k = 0; k < cnt; k++)
458 if (faceIds[n] == downFaces[k])
465 downFaces[cnt] = faceIds[n];
466 downTypes[cnt] = downTypesVol[n];
476 // ---------------------------------------------------------------------------
478 SMDS_Down2D::SMDS_Down2D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
479 SMDS_Downward(grid, nbDownCells)
482 _upCellTypes.clear();
487 SMDS_Down2D::~SMDS_Down2D()
491 int SMDS_Down2D::getNumberOfUpCells(int cellId)
494 if (_upCellIds[2 * cellId] >= 0)
496 if (_upCellIds[2 * cellId + 1] >= 0)
501 const int* SMDS_Down2D::getUpCells(int cellId)
503 //ASSERT((cellId >=0) && (cellId < _maxId));
504 return &_upCellIds[2 * cellId];
507 const unsigned char* SMDS_Down2D::getUpTypes(int cellId)
509 //ASSERT((cellId >=0) && (cellId < _maxId));
510 return &_upCellTypes[2 * cellId];
513 void SMDS_Down2D::getNodeIds(int cellId, std::set<int>& nodeSet)
515 for (int i = 0; i < _nbDownCells; i++)
517 int downCellId = _cellIds[_nbDownCells * cellId + i];
518 unsigned char cellType = _cellTypes[i];
519 this->_grid->getDownArray(cellType)->getNodeIds(downCellId, nodeSet);
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)
531 int SMDS_Down2D::computeVolumeIds(int cellId, int* ids)
533 // --- find point id's of the face
536 vtkIdType *pts; // will refer to the point id's of the face
537 _grid->GetCellPoints(cellId, npts, pts);
539 for (int i = 0; i < npts; i++)
540 nodes.push_back(pts[i]);
541 int nvol = this->computeVolumeIdsFromNodesFace(&nodes[0], npts, ids);
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.
550 * @param ids a couple of vtkId, initialized at -1 (no parent volume)
551 * @return number of volumes (0, 1 or 2)
553 int SMDS_Down2D::computeVolumeIds(ElemByNodesType& faceByNodes, int* ids)
555 int nvol = this->computeVolumeIdsFromNodesFace(&faceByNodes.nodeIds[0], faceByNodes.nbNodes, ids);
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
566 * @return number of volumes (0, 1 or 2)
568 int SMDS_Down2D::computeVolumeIdsFromNodesFace(int* pts, int npts, int* ids)
571 // --- find all the cells the points belong to, and how many of the points belong to a given cell
576 for (int i = 0; i < npts; i++)
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++)
584 int vtkCellId = cells[j];
586 for (int k = 0; k < cnt; k++)
588 if (cellIds[k] == vtkCellId)
597 cellIds[cnt] = vtkCellId;
599 // TODO ASSERT(cnt<1000);
605 // --- find the volume cells: they contains all the points and are of type volume
608 for (int i = 0; i < cnt; i++)
610 //MESSAGE("cell " << cellIds[i] << " points " << cellCnt[i]);
611 if (cellCnt[i] == npts)
613 int vtkElemId = cellIds[i];
614 int vtkType = _grid->GetCellType(vtkElemId);
615 if (SMDS_Downward::getCellDimension(vtkType) == 3)
617 ids[nvol] = vtkElemId; // store the volume id in given vector
628 void SMDS_Down2D::setTempNodes(int cellId, int vtkId)
631 vtkIdType *pts; // 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++)
637 _tempNodes[_nbNodes * cellId + i] = pts[i];
641 void SMDS_Down2D::setTempNodes(int cellId, ElemByNodesType& faceByNodes)
643 for (int i = 0; i < faceByNodes.nbNodes; i++)
644 _tempNodes[_nbNodes * cellId + i] = faceByNodes.nodeIds[i];
647 /*! Find if all the nodes belongs to the face.
649 * @param cellId the face cell Id
650 * @param nodeSet set of node id's to be found in the face list of nodes
653 bool SMDS_Down2D::isInFace(int cellId, int *pts, int npts)
656 int *nodes = &_tempNodes[_nbNodes * cellId];
657 for (int j = 0; j < npts; j++)
660 for (int i = 0; i < _nbNodes; i++)
662 if (nodes[i] == point)
669 return (nbFound == npts);
672 /*! Resize the downward connectivity storage vector if needed.
674 * @param nbElems total number of elements of the same type required
676 void SMDS_Down2D::allocate(int nbElems)
678 if (nbElems >= (int)_vtkCellIds.size())
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);
688 void SMDS_Down2D::compactStorage()
690 _cellIds.resize(_nbDownCells * _maxId);
691 _upCellIds.resize(2 * _maxId);
692 _upCellTypes.resize(2 * _maxId);
693 _vtkCellIds.resize(_maxId);
697 void SMDS_Down2D::addUpCell(int cellId, int upCellId, unsigned char aType)
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++)
706 vols[i] = upCellId; // use non affected volume
710 if ((vols[i] == upCellId) && (types[i] == aType)) // already done
716 int SMDS_Down2D::getNodeSet(int cellId, int* nodeSet)
718 for (int i = 0; i < _nbNodes; i++)
719 nodeSet[i] = _tempNodes[_nbNodes * cellId + i];
723 int SMDS_Down2D::FindEdgeByNodes(int cellId, ElemByNodesType& edgeByNodes)
725 int *edges = &_cellIds[_nbDownCells * cellId];
726 for (int i = 0; i < _nbDownCells; i++)
728 if ((edges[i] >= 0) && (edgeByNodes.vtkType == _cellTypes[i]))
731 int npts = this->_grid->getDownArray(edgeByNodes.vtkType)->getNodeSet(edges[i], nodeSet);
733 for (int j = 0; j < npts; j++)
735 int point = edgeByNodes.nodeIds[j];
737 for (int k = 0; k < npts; k++)
739 if (nodeSet[k] == point)
755 // ---------------------------------------------------------------------------
757 SMDS_Down3D::SMDS_Down3D(SMDS_UnstructuredGrid *grid, int nbDownCells) :
758 SMDS_Downward(grid, nbDownCells)
762 SMDS_Down3D::~SMDS_Down3D()
766 void SMDS_Down3D::allocate(int nbElems)
768 if (nbElems >= (int)_vtkCellIds.size())
770 _cellIds.resize(_nbDownCells * (nbElems + SMDS_Mesh::chunkSize), -1);
771 _vtkCellIds.resize(nbElems + SMDS_Mesh::chunkSize, -1);
775 void SMDS_Down3D::compactStorage()
777 // nothing to do, size was known before
780 int SMDS_Down3D::getNumberOfUpCells(int cellId)
785 const int* SMDS_Down3D::getUpCells(int cellId)
790 const unsigned char* SMDS_Down3D::getUpTypes(int cellId)
795 void SMDS_Down3D::getNodeIds(int cellId, std::set<int>& nodeSet)
797 int vtkId = this->_vtkCellIds[cellId];
799 vtkIdType *nodes; // 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]);
805 int SMDS_Down3D::FindFaceByNodes(int cellId, ElemByNodesType& faceByNodes)
807 int *faces = &_cellIds[_nbDownCells * cellId];
810 for (int i = 0; i < _nbDownCells; i++)
812 if ((faces[i] >= 0) && (faceByNodes.vtkType == _cellTypes[i]))
815 npoints = faceByNodes.nbNodes;
818 int npts = this->_grid->getDownArray(faceByNodes.vtkType)->getNodeSet(faces[i], nodeSet);
820 continue; // skip this face
822 for (int j = 0; j < npts; j++)
824 int point = faceByNodes.nodeIds[j];
826 for (int k = 0; k < npts; k++)
828 if (nodeSet[k] == point)
831 break; // point j is in the 2 faces, skip remaining k values
835 break; // point j is not in the 2 faces, skip the remaining tests
844 // ---------------------------------------------------------------------------
846 SMDS_DownEdge::SMDS_DownEdge(SMDS_UnstructuredGrid *grid) :
849 _cellTypes.push_back(VTK_VERTEX);
850 _cellTypes.push_back(VTK_VERTEX);
853 SMDS_DownEdge::~SMDS_DownEdge()
857 // ---------------------------------------------------------------------------
859 SMDS_DownQuadEdge::SMDS_DownQuadEdge(SMDS_UnstructuredGrid *grid) :
862 _cellTypes.push_back(VTK_VERTEX);
863 _cellTypes.push_back(VTK_VERTEX);
864 _cellTypes.push_back(VTK_VERTEX);
867 SMDS_DownQuadEdge::~SMDS_DownQuadEdge()
871 // ---------------------------------------------------------------------------
873 SMDS_DownTriangle::SMDS_DownTriangle(SMDS_UnstructuredGrid *grid) :
876 _cellTypes.push_back(VTK_LINE);
877 _cellTypes.push_back(VTK_LINE);
878 _cellTypes.push_back(VTK_LINE);
882 SMDS_DownTriangle::~SMDS_DownTriangle()
886 void SMDS_DownTriangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
888 int *nodes = &_tempNodes[_nbNodes * cellId];
889 edgesWithNodes.nbElems = 3;
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;
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;
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;
907 void SMDS_DownTriangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
916 faces[i] = lowCellId;
919 if (faces[i] == lowCellId)
925 // ---------------------------------------------------------------------------
927 SMDS_DownQuadTriangle::SMDS_DownQuadTriangle(SMDS_UnstructuredGrid *grid) :
930 _cellTypes.push_back(VTK_QUADRATIC_EDGE);
931 _cellTypes.push_back(VTK_QUADRATIC_EDGE);
932 _cellTypes.push_back(VTK_QUADRATIC_EDGE);
936 SMDS_DownQuadTriangle::~SMDS_DownQuadTriangle()
940 void SMDS_DownQuadTriangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
942 int *nodes = &_tempNodes[_nbNodes * cellId];
943 edgesWithNodes.nbElems = 3;
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;
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;
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;
964 void SMDS_DownQuadTriangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
973 edges[i] = lowCellId;
976 if (edges[i] == lowCellId)
982 // ---------------------------------------------------------------------------
984 SMDS_DownQuadrangle::SMDS_DownQuadrangle(SMDS_UnstructuredGrid *grid) :
987 _cellTypes.push_back(VTK_LINE);
988 _cellTypes.push_back(VTK_LINE);
989 _cellTypes.push_back(VTK_LINE);
990 _cellTypes.push_back(VTK_LINE);
994 SMDS_DownQuadrangle::~SMDS_DownQuadrangle()
998 void SMDS_DownQuadrangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
1000 int *nodes = &_tempNodes[_nbNodes * cellId];
1001 edgesWithNodes.nbElems = 4;
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;
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;
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;
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;
1024 void SMDS_DownQuadrangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1033 faces[i] = lowCellId;
1036 if (faces[i] == lowCellId)
1042 // ---------------------------------------------------------------------------
1044 SMDS_DownQuadQuadrangle::SMDS_DownQuadQuadrangle(SMDS_UnstructuredGrid *grid) :
1045 SMDS_Down2D(grid, 4)
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);
1054 SMDS_DownQuadQuadrangle::~SMDS_DownQuadQuadrangle()
1058 void SMDS_DownQuadQuadrangle::computeEdgesWithNodes(int cellId, ListElemByNodesType& edgesWithNodes)
1060 int *nodes = &_tempNodes[_nbNodes * cellId];
1061 edgesWithNodes.nbElems = 4;
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;
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;
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;
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;
1088 void SMDS_DownQuadQuadrangle::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1097 faces[i] = lowCellId;
1100 if (faces[i] == lowCellId)
1106 // ---------------------------------------------------------------------------
1108 SMDS_DownTetra::SMDS_DownTetra(SMDS_UnstructuredGrid *grid) :
1109 SMDS_Down3D(grid, 4)
1111 _cellTypes.push_back(VTK_TRIANGLE);
1112 _cellTypes.push_back(VTK_TRIANGLE);
1113 _cellTypes.push_back(VTK_TRIANGLE);
1114 _cellTypes.push_back(VTK_TRIANGLE);
1117 SMDS_DownTetra::~SMDS_DownTetra()
1121 void SMDS_DownTetra::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1125 for ( size_t i = 0; i < orderedNodes.size(); i++ )
1126 setNodes.insert(orderedNodes[i]);
1127 //MESSAGE("cellId = " << cellId);
1130 vtkIdType *nodes; // will refer to the point id's of the volume
1131 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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++)
1139 for (int i = 0; i < 3; i++)
1140 tofind.insert(nodes[ids[3 * k + i]]);
1141 if (setNodes == tofind)
1143 for (int i = 0; i < 3; i++)
1144 orderedNodes[i] = nodes[ids[3 * k + i]];
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]);
1153 void SMDS_DownTetra::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1162 faces[i] = lowCellId;
1165 if (faces[i] == lowCellId)
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
1177 void SMDS_DownTetra::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1179 // --- find point id's of the volume
1182 vtkIdType *nodes; // will refer to the point id's of the volume
1183 _grid->GetCellPoints(cellId, npts, nodes);
1185 // --- create all the ordered list of node id's for each face
1187 facesWithNodes.nbElems = 4;
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;
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;
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;
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;
1214 // ---------------------------------------------------------------------------
1216 SMDS_DownQuadTetra::SMDS_DownQuadTetra(SMDS_UnstructuredGrid *grid) :
1217 SMDS_Down3D(grid, 4)
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);
1225 SMDS_DownQuadTetra::~SMDS_DownQuadTetra()
1229 void SMDS_DownQuadTetra::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1233 for ( size_t i = 0; i < orderedNodes.size(); i++ )
1234 setNodes.insert(orderedNodes[i]);
1235 //MESSAGE("cellId = " << cellId);
1238 vtkIdType *nodes; // will refer to the point id's of the volume
1239 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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++)
1247 for (int i = 0; i < 6; i++)
1248 tofind.insert(nodes[ids[6 * k + i]]);
1249 if (setNodes == tofind)
1251 for (int i = 0; i < 6; i++)
1252 orderedNodes[i] = nodes[ids[6 * k + i]];
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]);
1261 void SMDS_DownQuadTetra::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1270 faces[i] = lowCellId;
1273 if (faces[i] == lowCellId)
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
1287 void SMDS_DownQuadTetra::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1289 // --- find point id's of the volume
1292 vtkIdType *nodes; // will refer to the point id's of the volume
1293 _grid->GetCellPoints(cellId, npts, nodes);
1295 // --- create all the ordered list of node id's for each face
1297 facesWithNodes.nbElems = 4;
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;
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;
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;
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;
1336 // ---------------------------------------------------------------------------
1338 SMDS_DownPyramid::SMDS_DownPyramid(SMDS_UnstructuredGrid *grid) :
1339 SMDS_Down3D(grid, 5)
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);
1348 SMDS_DownPyramid::~SMDS_DownPyramid()
1352 void SMDS_DownPyramid::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1356 for ( size_t i = 0; i < orderedNodes.size(); i++ )
1357 setNodes.insert(orderedNodes[i]);
1358 //MESSAGE("cellId = " << cellId);
1361 vtkIdType *nodes; // will refer to the point id's of the volume
1362 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
1365 int ids[16] = { 0, 1, 2, 3, 0, 3, 4, 3, 2, 4, 2, 1, 4, 1, 0, 4 };
1367 // Quadrangular face
1369 for (int i = 0; i < 4; i++)
1370 tofind.insert(nodes[ids[i]]);
1371 if (setNodes == tofind)
1373 for (int i = 0; i < 4; i++)
1374 orderedNodes[i] = nodes[ids[i]];
1378 for (int k = 0; k < 4; k++)
1381 for (int i = 0; i < 3; i++)
1382 tofind.insert(nodes[ids[4 + 3 * k + i]]);
1383 if (setNodes == tofind)
1385 for (int i = 0; i < 3; i++)
1386 orderedNodes[i] = nodes[ids[4 + 3 * k + i]];
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]);
1395 void SMDS_DownPyramid::addDownCell(int cellId, int lowCellId, unsigned char aType)
1397 //ASSERT((cellId >=0) && (cellId < _maxId));
1398 int *faces = &_cellIds[_nbDownCells * cellId];
1399 if (aType == VTK_QUAD)
1403 faces[0] = lowCellId;
1406 if (faces[0] == lowCellId)
1411 //ASSERT(aType == VTK_TRIANGLE);
1412 for (int i = 1; i < _nbDownCells; i++)
1416 faces[i] = lowCellId;
1419 if (faces[i] == lowCellId)
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
1434 void SMDS_DownPyramid::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1436 // --- find point id's of the volume
1439 vtkIdType *nodes; // will refer to the point id's of the volume
1440 _grid->GetCellPoints(cellId, npts, nodes);
1442 // --- create all the ordered list of node id's for each face
1444 facesWithNodes.nbElems = 5;
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;
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;
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;
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;
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;
1478 // ---------------------------------------------------------------------------
1480 SMDS_DownQuadPyramid::SMDS_DownQuadPyramid(SMDS_UnstructuredGrid *grid) :
1481 SMDS_Down3D(grid, 5)
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);
1490 SMDS_DownQuadPyramid::~SMDS_DownQuadPyramid()
1494 void SMDS_DownQuadPyramid::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1496 // MESSAGE("SMDS_DownQuadPyramid::getOrderedNodesOfFace cellId = " << cellId);
1499 for ( size_t i = 0; i < orderedNodes.size(); i++)
1500 setNodes.insert(orderedNodes[i]);
1501 //MESSAGE("cellId = " << cellId);
1504 vtkIdType *nodes; // will refer to the point id's of the volume
1505 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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 };
1511 // Quadrangular face
1513 for (int i = 0; i < 8; i++)
1514 tofind.insert(nodes[ids[i]]);
1515 if (setNodes == tofind)
1517 for (int i = 0; i < 8; i++)
1518 orderedNodes[i] = nodes[ids[i]];
1522 for (int k = 0; k < 4; k++)
1525 for (int i = 0; i < 6; i++)
1526 tofind.insert(nodes[ids[8 + 6 * k + i]]);
1527 if (setNodes == tofind)
1529 for (int i = 0; i < 6; i++)
1530 orderedNodes[i] = nodes[ids[8 + 6 * k + i]];
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]);
1539 void SMDS_DownQuadPyramid::addDownCell(int cellId, int lowCellId, unsigned char aType)
1541 //ASSERT((cellId >=0) && (cellId < _maxId));
1542 int *faces = &_cellIds[_nbDownCells * cellId];
1543 if (aType == VTK_QUADRATIC_QUAD)
1547 faces[0] = lowCellId;
1550 if (faces[0] == lowCellId)
1555 //ASSERT(aType == VTK_QUADRATIC_TRIANGLE);
1556 for (int i = 1; i < _nbDownCells; i++)
1560 faces[i] = lowCellId;
1563 if (faces[i] == lowCellId)
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
1579 void SMDS_DownQuadPyramid::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1581 // --- find point id's of the volume
1584 vtkIdType *nodes; // will refer to the point id's of the volume
1585 _grid->GetCellPoints(cellId, npts, nodes);
1587 // --- create all the ordered list of node id's for each face
1589 facesWithNodes.nbElems = 5;
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;
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;
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;
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;
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;
1639 // ---------------------------------------------------------------------------
1641 SMDS_DownPenta::SMDS_DownPenta(SMDS_UnstructuredGrid *grid) :
1642 SMDS_Down3D(grid, 5)
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);
1651 SMDS_DownPenta::~SMDS_DownPenta()
1655 void SMDS_DownPenta::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1659 for ( size_t i = 0; i < orderedNodes.size(); i++)
1660 setNodes.insert(orderedNodes[i]);
1661 //MESSAGE("cellId = " << cellId);
1664 vtkIdType *nodes; // will refer to the point id's of the volume
1665 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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 };
1672 for (int k = 0; k < 2; k++)
1675 for (int i = 0; i < 3; i++)
1676 tofind.insert(nodes[ids[3 * k + i]]);
1677 if (setNodes == tofind)
1679 for (int i = 0; i < 3; i++)
1680 orderedNodes[i] = nodes[ids[3 * k + i]];
1684 // Quadrangular faces
1685 for (int k = 0; k < 3; k++)
1688 for (int i = 0; i < 4; i++)
1689 tofind.insert(nodes[ids[6 + 4 * k + i]]);
1690 if (setNodes == tofind)
1692 for (int i = 0; i < 4; i++)
1693 orderedNodes[i] = nodes[ids[6 + 4 * k + i]];
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]);
1702 void SMDS_DownPenta::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1711 faces[i] = lowCellId;
1714 if (faces[i] == lowCellId)
1719 //ASSERT(aType == VTK_TRIANGLE);
1720 for (int i = 3; i < _nbDownCells; i++)
1724 faces[i] = lowCellId;
1727 if (faces[i] == lowCellId)
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
1743 void SMDS_DownPenta::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1745 // --- find point id's of the volume
1748 vtkIdType *nodes; // will refer to the point id's of the volume
1749 _grid->GetCellPoints(cellId, npts, nodes);
1751 // --- create all the ordered list of node id's for each face
1753 facesWithNodes.nbElems = 5;
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;
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;
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;
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;
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;
1789 // ---------------------------------------------------------------------------
1791 SMDS_DownQuadPenta::SMDS_DownQuadPenta(SMDS_UnstructuredGrid *grid) :
1792 SMDS_Down3D(grid, 5)
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);
1801 SMDS_DownQuadPenta::~SMDS_DownQuadPenta()
1805 void SMDS_DownQuadPenta::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1809 for ( size_t i = 0; i < orderedNodes.size(); i++ )
1810 setNodes.insert(orderedNodes[i]);
1811 //MESSAGE("cellId = " << cellId);
1814 vtkIdType *nodes; // will refer to the point id's of the volume
1815 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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 };
1823 for (int k = 0; k < 2; k++)
1826 for (int i = 0; i < 6; i++)
1827 tofind.insert(nodes[ids[6 * k + i]]);
1828 if (setNodes == tofind)
1830 for (int i = 0; i < 6; i++)
1831 orderedNodes[i] = nodes[ids[6 * k + i]];
1835 // Quadrangular faces
1836 for (int k = 0; k < 3; k++)
1839 for (int i = 0; i < 8; i++)
1840 tofind.insert(nodes[ids[12 + 8 * k + i]]);
1841 if (setNodes == tofind)
1843 for (int i = 0; i < 8; i++)
1844 orderedNodes[i] = nodes[ids[12 + 8 * k + i]];
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]);
1853 void SMDS_DownQuadPenta::addDownCell(int cellId, int lowCellId, unsigned char aType)
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++)
1862 faces[i] = lowCellId;
1865 if (faces[i] == lowCellId)
1870 //ASSERT(aType == VTK_QUADRATIC_TRIANGLE);
1871 for (int i = 3; i < _nbDownCells; i++)
1875 faces[i] = lowCellId;
1878 if (faces[i] == lowCellId)
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
1895 void SMDS_DownQuadPenta::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
1897 // --- find point id's of the volume
1900 vtkIdType *nodes; // will refer to the point id's of the volume
1901 _grid->GetCellPoints(cellId, npts, nodes);
1903 // --- create all the ordered list of node id's for each face
1905 facesWithNodes.nbElems = 5;
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;
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;
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;
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;
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;
1959 // ---------------------------------------------------------------------------
1961 SMDS_DownHexa::SMDS_DownHexa(SMDS_UnstructuredGrid *grid) :
1962 SMDS_Down3D(grid, 6)
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);
1972 SMDS_DownHexa::~SMDS_DownHexa()
1976 void SMDS_DownHexa::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
1980 for ( size_t i = 0; i < orderedNodes.size(); i++ )
1981 setNodes.insert(orderedNodes[i]);
1982 //MESSAGE("cellId = " << cellId);
1985 vtkIdType *nodes; // will refer to the point id's of the volume
1986 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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
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)
1998 for (int i = 0; i < 4; i++)
1999 orderedNodes[i] = nodes[ids[4 * k + i]];
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]);
2009 void SMDS_DownHexa::addDownCell(int cellId, int lowCellId, unsigned char aType)
2011 //ASSERT((cellId >=0)&& (cellId < _maxId));
2012 int *faces = &_cellIds[_nbDownCells * cellId];
2013 for (int i = 0; i < _nbDownCells; i++)
2017 faces[i] = lowCellId;
2020 if (faces[i] == lowCellId)
2024 // MESSAGE("-------------------------------------> trop de faces ! " << cellId << " " << lowCellId);
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
2035 void SMDS_DownHexa::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
2037 // --- find point id's of the volume
2040 vtkIdType *nodes; // will refer to the point id's of the volume
2041 _grid->GetCellPoints(cellId, npts, nodes);
2043 // --- create all the ordered list of node id's for each face
2045 facesWithNodes.nbElems = 6;
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;
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;
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;
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;
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;
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;
2090 // ---------------------------------------------------------------------------
2092 SMDS_DownQuadHexa::SMDS_DownQuadHexa(SMDS_UnstructuredGrid *grid) :
2093 SMDS_Down3D(grid, 6)
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);
2103 SMDS_DownQuadHexa::~SMDS_DownQuadHexa()
2107 void SMDS_DownQuadHexa::getOrderedNodesOfFace(int cellId, std::vector<vtkIdType>& orderedNodes)
2111 for ( size_t i = 0; i < orderedNodes.size(); i++ )
2112 setNodes.insert(orderedNodes[i]);
2113 //MESSAGE("cellId = " << cellId);
2116 vtkIdType *nodes; // will refer to the point id's of the volume
2117 _grid->GetCellPoints(this->_vtkCellIds[cellId], npts, nodes);
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++)
2126 for (int i = 0; i < 8; i++)
2127 tofind.insert(nodes[ids[8 * k + i]]);
2128 if (setNodes == tofind)
2130 for (int i = 0; i < 8; i++)
2131 orderedNodes[i] = nodes[ids[8 * k + i]];
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]);
2140 void SMDS_DownQuadHexa::addDownCell(int cellId, int lowCellId, unsigned char aType)
2142 //ASSERT((cellId >=0)&& (cellId < _maxId));
2143 int *faces = &_cellIds[_nbDownCells * cellId];
2144 for (int i = 0; i < _nbDownCells; i++)
2148 faces[i] = lowCellId;
2151 if (faces[i] == lowCellId)
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
2166 void SMDS_DownQuadHexa::computeFacesWithNodes(int cellId, ListElemByNodesType& facesWithNodes)
2168 // --- find point id's of the volume
2171 vtkIdType *nodes; // will refer to the point id's of the volume
2172 _grid->GetCellPoints(cellId, npts, nodes);
2174 // --- create all the ordered list of node id's for each face
2176 facesWithNodes.nbElems = 6;
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;
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;
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;
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;
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;
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;
2245 // ---------------------------------------------------------------------------