Chi-Tech
unpartmesh_00c_vtk_cellsvtk2chi.cc
Go to the documentation of this file.
2
3#include <vtkPolygon.h>
4#include <vtkLine.h>
5#include <vtkVertex.h>
6
7// ###################################################################
8/**Creates a raw polyhedron cell from a vtk-polyhedron.*/
11{
12 const std::string fname =
13 "chi_mesh::UnpartitionedMesh::CreateCellFromVTKPolyhedron";
14
15 CellType sub_type;
16 switch (vtk_cell->GetCellType())
17 {
18 case VTK_HEXAGONAL_PRISM:
19 case VTK_PENTAGONAL_PRISM:
20 case VTK_POLYHEDRON:
21 sub_type = CellType::POLYHEDRON;
22 break;
23 case VTK_PYRAMID:
24 sub_type = CellType::PYRAMID;
25 break;
26 case VTK_WEDGE:
27 sub_type = CellType::WEDGE;
28 break;
29 case VTK_HEXAHEDRON:
30 case VTK_VOXEL:
31 sub_type = CellType::HEXAHEDRON;
32 break;
33 case VTK_TETRA:
34 sub_type = CellType::TETRAHEDRON;
35 break;
36 default:
37 throw std::logic_error(fname + ": Unsupported 3D cell type encountered.");
38 }
39 auto polyh_cell = new LightWeightCell(CellType::POLYHEDRON, sub_type);
40
41 auto num_cpoints = vtk_cell->GetNumberOfPoints();
42 auto num_cfaces = vtk_cell->GetNumberOfFaces();
43
44 polyh_cell->vertex_ids.reserve(num_cpoints);
45 auto point_ids = vtk_cell->GetPointIds();
46 for (int p = 0; p < num_cpoints; ++p)
47 {
48 uint64_t point_id = point_ids->GetId(p);
49 polyh_cell->vertex_ids.push_back(point_id);
50 } // for p
51
52 switch (sub_type)
53 {
54 // The cell vertex ids in VTK is the same as in ChiTech so we don't
55 // need to remap the vertices. We do however need to remap the faces.
57 {
58 std::vector<std::vector<uint64_t>> face_vids = {{1, 2, 6, 5},
59 {3, 0, 4, 7},
60 {2, 3, 7, 6},
61 {0, 1, 5, 4},
62 {4, 5, 6, 7},
63 {3, 2, 1, 0}};
64 for (int f = 0; f < 6; ++f)
65 {
66 LightWeightFace face;
67
68 face.vertex_ids.reserve(4);
69 for (int p = 0; p < 4; ++p)
70 face.vertex_ids.push_back(polyh_cell->vertex_ids[face_vids[f][p]]);
71
72 polyh_cell->faces.push_back(face);
73 }
74 break;
75 }
76 // For wedges we need to remap cell vertices and faces.
77 case CellType::WEDGE:
78 {
79 std::vector<uint64_t> remapping = {0, 2, 1, 3, 5, 4};
80 std::vector<uint64_t> remapped(6, 0);
81 for (int i = 0; i < 6; ++i)
82 remapped[i] = polyh_cell->vertex_ids[remapping[i]];
83 polyh_cell->vertex_ids = remapped;
84
85 std::vector<std::vector<uint64_t>> face_vids = {
86 {0, 1, 4, 3}, {1, 2, 5, 4}, {2, 0, 3, 5}, {3, 4, 5}, {0, 2, 1}};
87 for (int f = 0; f < 6; ++f)
88 {
89 LightWeightFace face;
90
91 face.vertex_ids.reserve(4);
92 for (int p = 0; p < face_vids[f].size(); ++p)
93 face.vertex_ids.push_back(polyh_cell->vertex_ids[face_vids[f][p]]);
94
95 polyh_cell->faces.push_back(face);
96 }
97 break;
98 }
99 // We dont need cell vertex id remapping but we need to impose our own
100 // face orientation
102 {
103 std::vector<std::vector<uint64_t>> face_vids = {{0, 2, 1},
104 {0, 1, 3},
105 {0, 3, 2},
106 {3, 1, 2}};
107 for (int f = 0; f < 4; ++f)
108 {
109 LightWeightFace face;
110
111 face.vertex_ids.reserve(3);
112 for (int p = 0; p < 3; ++p)
113 face.vertex_ids.push_back(polyh_cell->vertex_ids[face_vids[f][p]]);
114
115 polyh_cell->faces.push_back(face);
116 }
117 break;
118 }
119 default:
120 {
121 polyh_cell->faces.reserve(num_cfaces);
122 for (int f = 0; f < num_cfaces; ++f)
123 {
124 LightWeightFace face;
125 auto vtk_face = vtk_cell->GetFace(f);
126 auto num_face_points = vtk_face->GetNumberOfPoints();
127
128 face.vertex_ids.reserve(num_face_points);
129 auto face_point_ids = vtk_face->GetPointIds();
130 for (int p = 0; p < num_face_points; ++p)
131 {
132 uint64_t point_id = face_point_ids->GetId(p);
133 face.vertex_ids.push_back(point_id);
134 }
135
136 polyh_cell->faces.push_back(face);
137 }
138 break;
139 } // default
140 } // switch on sub_type
141
142 return polyh_cell;
143}
144
145// ###################################################################
146/**Creates a raw polygon cell from a vtk-polygon.*/
149{
150 const std::string fname =
151 "chi_mesh::UnpartitionedMesh::CreateCellFromVTKPolygon";
152
153 CellType sub_type;
154 switch (vtk_cell->GetCellType())
155 {
156 case VTK_POLYGON:
157 sub_type = CellType::POLYGON;
158 break;
159 case VTK_QUAD:
160 case VTK_PIXEL:
161 sub_type = CellType::QUADRILATERAL;
162 break;
163 case VTK_TRIANGLE:
164 sub_type = CellType::TRIANGLE;
165 break;
166 default:
167 throw std::logic_error(fname + ": Unsupported 2D cell type encountered.");
168 }
169
170 auto poly_cell = new LightWeightCell(CellType::POLYGON, sub_type);
171
172 auto num_cpoints = vtk_cell->GetNumberOfPoints();
173 auto num_cfaces = num_cpoints;
174
175 poly_cell->vertex_ids.reserve(num_cpoints);
176 auto point_ids = vtk_cell->GetPointIds();
177 for (int p = 0; p < num_cpoints; ++p)
178 {
179 uint64_t point_id = point_ids->GetId(p);
180 poly_cell->vertex_ids.push_back(point_id);
181 } // for p
182
183 poly_cell->faces.reserve(num_cfaces);
184 for (int f = 0; f < num_cfaces; ++f)
185 {
186 LightWeightFace face;
187
188 auto v0_id = poly_cell->vertex_ids[f];
189 auto v1_id = (f < (num_cfaces - 1)) ? poly_cell->vertex_ids[f + 1]
190 : poly_cell->vertex_ids[0];
191
192 face.vertex_ids.reserve(2);
193 face.vertex_ids.push_back(v0_id);
194 face.vertex_ids.push_back(v1_id);
195
196 poly_cell->faces.push_back(face);
197 }
198
199 return poly_cell;
200}
201
202// ###################################################################
203/**Creates a raw slab cell from a vtk-line.*/
206{
207 const std::string fname =
208 "chi_mesh::UnpartitionedMesh::CreateCellFromVTKPolygon";
209
210 CellType sub_type;
211 switch (vtk_cell->GetCellType())
212 {
213 case VTK_LINE:
214 sub_type = CellType::SLAB;
215 break;
216 default:
217 throw std::logic_error(fname + ": Unsupported 1D cell type encountered.");
218 }
219
220 auto slab_cell = new LightWeightCell(CellType::SLAB, sub_type);
221
222 auto vtk_line = vtkLine::SafeDownCast(vtk_cell);
223 auto num_cpoints = vtk_line->GetNumberOfPoints();
224 auto num_cfaces = num_cpoints;
225
226 slab_cell->vertex_ids.reserve(num_cpoints);
227 auto point_ids = vtk_line->GetPointIds();
228 for (int p = 0; p < num_cpoints; ++p)
229 {
230 uint64_t point_id = point_ids->GetId(p);
231 slab_cell->vertex_ids.push_back(point_id);
232 } // for p
233
234 slab_cell->faces.reserve(num_cfaces);
235 for (int f = 0; f < num_cfaces; ++f)
236 {
237 LightWeightFace face;
238
239 auto v_id = slab_cell->vertex_ids[f];
240
241 face.vertex_ids.reserve(1);
242 face.vertex_ids.push_back(v_id);
243
244 slab_cell->faces.push_back(face);
245 }
246
247 return slab_cell;
248}
249
250// ###################################################################
251/**Creates a raw point cell from a vtk-vertex.*/
254{
255 auto point_cell = new LightWeightCell(CellType::GHOST, CellType::POINT);
256
257 auto vtk_vertex = vtkVertex::SafeDownCast(vtk_cell);
258 auto num_cpoints = vtk_vertex->GetNumberOfPoints();
259
260 point_cell->vertex_ids.reserve(num_cpoints);
261 auto point_ids = vtk_vertex->GetPointIds();
262 for (int p = 0; p < num_cpoints; ++p)
263 {
264 uint64_t point_id = point_ids->GetId(p);
265 point_cell->vertex_ids.push_back(point_id);
266 } // for p
267
268 return point_cell;
269}
static LightWeightCell * CreateCellFromVTKLine(vtkCell *vtk_cell)
static LightWeightCell * CreateCellFromVTKVertex(vtkCell *vtk_cell)
static LightWeightCell * CreateCellFromVTKPolygon(vtkCell *vtk_cell)
static LightWeightCell * CreateCellFromVTKPolyhedron(vtkCell *vtk_cell)
CellType
Definition: cell.h:12