Salome HOME
bos #20256: [CEA 18523] Porting SMESH to int 64 bits
[modules/smesh.git] / src / SMDS / SMDS_VtkCellIterator.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 #include "SMDS_VtkCellIterator.hxx"
21 #include "utilities.h"
22
23 #include <vtkCell.h>
24 #include <vtkIdList.h>
25
26 _GetVtkNodes::_GetVtkNodes( TVtkIdList&        vtkIds,
27                             SMDS_Mesh*         mesh,
28                             vtkIdType          vtkCellId,
29                             SMDSAbs_EntityType type )
30 {
31   vtkUnstructuredGrid*         grid = mesh->GetGrid();
32   const std::vector<int>& interlace = SMDS_MeshCell::fromVtkOrder( type );
33   vtkIdType npts;
34   vtkIdType const *pts(nullptr);
35   grid->GetCellPoints( vtkCellId, npts, pts );
36   vtkIds.resize( npts );
37   if ( interlace.empty() )
38   {
39     vtkIds.assign( pts, pts + npts );
40   }
41   else
42   {
43     for (vtkIdType i = 0; i < npts; i++)
44       vtkIds[ i ] = pts[ interlace[i] ];
45   }
46 }
47
48 _GetVtkNodesToUNV::_GetVtkNodesToUNV( TVtkIdList&        vtkIds,
49                                       SMDS_Mesh*         mesh,
50                                       vtkIdType          vtkCellId,
51                                       SMDSAbs_EntityType type )
52 {
53   vtkUnstructuredGrid* grid = mesh->GetGrid();
54   vtkIdType npts;
55   vtkIdType const *pts(nullptr);
56   grid->GetCellPoints( vtkCellId, npts, pts );
57   const int *ids = 0;
58   switch ( type )
59   {
60   case SMDSEntity_Quad_Edge:
61   {
62     static int id[] = { 0, 2, 1 };
63     ids = id;
64     break;
65   }
66   case SMDSEntity_Quad_Triangle:
67   case SMDSEntity_BiQuad_Triangle:
68   {
69     static int id[] = { 0, 3, 1, 4, 2, 5 };
70     ids = id;
71     npts = 6;
72     break;
73   }
74   case SMDSEntity_Quad_Quadrangle:
75   case SMDSEntity_BiQuad_Quadrangle:
76   {
77     static int id[] = { 0, 4, 1, 5, 2, 6, 3, 7 };
78     ids = id;
79     npts = 8;
80     break;
81   }
82   case SMDSEntity_Quad_Tetra:
83   {
84     static int id[] = { 0, 4, 1, 5, 2, 6, 7, 8, 9, 3 };
85     ids = id;
86     break;
87   }
88   case SMDSEntity_Quad_Pyramid:
89   {
90     static int id[] = { 0, 5, 1, 6, 2, 7, 3, 8, 9, 10, 11, 12, 4 };
91     ids = id;
92     break;
93   }
94   case SMDSEntity_Penta:
95   {
96     static int id[] = { 0, 2, 1, 3, 5, 4 };
97     ids = id;
98     break;
99   }
100   case SMDSEntity_Quad_Penta:
101   case SMDSEntity_BiQuad_Penta: //TODO: check
102   {
103     static int id[] = { 0, 8, 2, 7, 1, 6, 12, 14, 13, 3, 11, 5, 10, 4, 9 };
104     ids = id;
105     break;
106   }
107   case SMDSEntity_Quad_Hexa:
108   case SMDSEntity_TriQuad_Hexa:
109   {
110     static int id[] = { 0, 8, 1, 9, 2, 10, 3, 11, 16, 17, 18, 19, 4, 12, 5, 13, 6, 14, 7, 15 };
111     ids = id;
112     npts = 20;
113     break;
114   }
115   case SMDSEntity_Polygon:
116   case SMDSEntity_Quad_Polygon:
117   case SMDSEntity_Polyhedra:
118   case SMDSEntity_Quad_Polyhedra:
119   default:
120     const std::vector<int>& i = SMDS_MeshCell::interlacedSmdsOrder( type, npts );
121     if ( !i.empty() )
122       ids = & i[0];
123   }
124
125   vtkIds.resize( npts );
126
127   if ( ids )
128     for (int i = 0; i < npts; i++)
129       vtkIds[ i ] =  pts[ ids[i] ];
130   else
131     vtkIds.assign( pts, pts + npts );
132 }
133
134 _GetVtkNodesPolyh::_GetVtkNodesPolyh( TVtkIdList&        vtkIds,
135                                       SMDS_Mesh*         mesh,
136                                       vtkIdType          vtkCellId,
137                                       SMDSAbs_EntityType type )
138 {
139   vtkUnstructuredGrid* grid = mesh->GetGrid();
140   switch ( type )
141   {
142   case SMDSEntity_Polyhedra:
143   {
144     vtkIdType nFaces = 0;
145     vtkIdType const *ptIds(nullptr);
146     grid->GetFaceStream( vtkCellId, nFaces, ptIds );
147     int id = 0, nbNodesInFaces = 0;
148     for ( int i = 0; i < nFaces; i++ )
149     {
150       int nodesInFace = ptIds[id]; // nodeIds in ptIds[id+1 .. id+nodesInFace]
151       nbNodesInFaces += nodesInFace;
152       id += (nodesInFace + 1);
153     }
154     vtkIds.resize( nbNodesInFaces );
155     id = 0;
156     int n = 0;
157     for ( int i = 0; i < nFaces; i++ )
158     {
159       int nodesInFace = ptIds[id]; // nodeIds in ptIds[id+1 .. id+nodesInFace]
160       for ( int k = 1; k <= nodesInFace; k++ )
161         vtkIds[ n++ ] = ptIds[ id + k ];
162       id += (nodesInFace + 1);
163     }
164     break;
165   }
166   default:
167     assert(0);
168   }
169 }