Chi-Tech
chi_grid_vtk_utils_00_celluploading.cc
Go to the documentation of this file.
2
3#include "chi_meshcontinuum.h"
4
5#include <vtkPoints.h>
6#include <vtkUnstructuredGrid.h>
7
8//###################################################################
9/**Uploads vertices and cells to an unstructured grid.*/
11 const chi_mesh::Cell &cell,
12 int64_t &node_counter,
13 vtkNew<vtkPoints> &points,
15{
16 size_t num_verts = cell.vertex_ids_.size();
17
18 std::vector<vtkIdType> cell_vids(num_verts);
19 for (size_t v=0; v<num_verts; v++)
20 {
21 uint64_t vgi = cell.vertex_ids_[v];
22 std::vector<double> d_node(3);
23 d_node[0] = grid.vertices[vgi].x;
24 d_node[1] = grid.vertices[vgi].y;
25 d_node[2] = grid.vertices[vgi].z;
26
27 points->InsertPoint(node_counter,d_node.data());
28 cell_vids[v] = node_counter++;
29 }
30
31 if (cell.Type() == chi_mesh::CellType::SLAB)
32 {
33 ugrid->InsertNextCell(VTK_LINE,
34 static_cast<vtkIdType>(num_verts),
35 cell_vids.data());
36 }
38 {
39 int vtk_subtype;
40 switch (cell.SubType())
41 {
42 case CellType::POLYGON: vtk_subtype = VTK_POLYGON; break;
43 case CellType::QUADRILATERAL: vtk_subtype = VTK_QUAD; break;
44 case CellType::TRIANGLE: vtk_subtype = VTK_TRIANGLE; break;
45 default: vtk_subtype = VTK_POLYGON; break;
46 }
47
48 ugrid->InsertNextCell(vtk_subtype,
49 static_cast<vtkIdType>(num_verts),
50 cell_vids.data());
51 }
53 {
54 // Build polyhedron faces
55 std::vector<vtkIdType> faces_vids;
56
57 size_t num_faces = cell.faces_.size();
58 for (auto& face : cell.faces_)
59 {
60 size_t num_fverts = face.vertex_ids_.size();
61 std::vector<vtkIdType> face_info(num_fverts);
62 for (size_t fv=0; fv<num_fverts; fv++)
63 {
64 size_t v = 0;
65 for (size_t cv=0; cv<num_verts; ++cv)
66 if (cell.vertex_ids_[cv] == face.vertex_ids_[fv])
67 { v = cv; break; }
68
69 face_info[fv] = cell_vids[v];
70 }
71
72 faces_vids.push_back(static_cast<vtkIdType>(num_fverts));
73 for (auto vid : face_info)
74 faces_vids.push_back(vid);
75 }//for f
76
77 int vtk_subtype;
78 switch (cell.SubType())
79 {
80 case CellType::POLYHEDRON: vtk_subtype = VTK_POLYHEDRON; break;
81 case CellType::PYRAMID: vtk_subtype = VTK_PYRAMID; break;
82 case CellType::WEDGE: vtk_subtype = VTK_WEDGE; break;
83 case CellType::HEXAHEDRON: vtk_subtype = VTK_HEXAHEDRON; break;
84 case CellType::TETRAHEDRON: vtk_subtype = VTK_TETRA; break;
85 default: vtk_subtype = VTK_POLYHEDRON; break;
86 }
87
88 ugrid->InsertNextCell(vtk_subtype,
89 static_cast<vtkIdType>(num_verts),
90 cell_vids.data(),
91 static_cast<vtkIdType>(num_faces),
92 faces_vids.data());
93 }//polyhedron
94}
95
96//###################################################################
97/**Uploads vertices and cells to an unstructured grid.*/
100 const std::vector<uint64_t>& vertex_map,
102{
103 size_t num_verts = cell.vertex_ids_.size();
104
105 std::vector<vtkIdType> cell_vids(num_verts);
106 for (size_t v=0; v<num_verts; v++)
107 cell_vids[v] = static_cast<vtkIdType>(vertex_map[cell.vertex_ids_[v]]);
108
109 if (cell.Type() == chi_mesh::CellType::SLAB)
110 {
111 ugrid->InsertNextCell(VTK_LINE,
112 static_cast<vtkIdType>(num_verts),
113 cell_vids.data());
114 }
115 if (cell.Type() == chi_mesh::CellType::POLYGON)
116 {
117 int vtk_subtype;
118 switch (cell.SubType())
119 {
120 case CellType::POLYGON: vtk_subtype = VTK_POLYGON; break;
121 case CellType::QUADRILATERAL: vtk_subtype = VTK_QUAD; break;
122 case CellType::TRIANGLE: vtk_subtype = VTK_TRIANGLE; break;
123 default: vtk_subtype = VTK_POLYGON; break;
124 }
125
126 ugrid->InsertNextCell(vtk_subtype,
127 static_cast<vtkIdType>(num_verts),
128 cell_vids.data());
129 }
131 {
132 // Build polyhedron faces
133 std::vector<vtkIdType> faces_vids;
134
135 int vtk_subtype;
136 switch (cell.SubType())
137 {
138 case CellType::POLYHEDRON: vtk_subtype = VTK_POLYHEDRON; break;
139 case CellType::PYRAMID: vtk_subtype = VTK_PYRAMID; break;
140 case CellType::WEDGE: vtk_subtype = VTK_WEDGE; break;
141 case CellType::HEXAHEDRON: vtk_subtype = VTK_HEXAHEDRON; break;
142 case CellType::TETRAHEDRON: vtk_subtype = VTK_TETRA; break;
143 default: vtk_subtype = VTK_POLYHEDRON; break;
144 }
145
146 switch (cell.SubType())
147 {
148 case CellType::POLYHEDRON:
149 {
150 size_t num_faces = cell.faces_.size();
151 for (auto &face: cell.faces_)
152 {
153 size_t num_fverts = face.vertex_ids_.size();
154 std::vector<vtkIdType> face_info(num_fverts);
155 for (size_t fv = 0; fv < num_fverts; fv++)
156 {
157 size_t v = 0;
158 for (size_t cv = 0; cv < num_verts; ++cv)
159 if (cell.vertex_ids_[cv] == face.vertex_ids_[fv])
160 {
161 v = cv;
162 break;
163 }
164
165 face_info[fv] = cell_vids[v];
166 }
167
168 faces_vids.push_back(static_cast<vtkIdType>(num_fverts));
169 for (auto vid: face_info)
170 faces_vids.push_back(vid);
171 }//for f
172
173 ugrid->InsertNextCell(vtk_subtype,
174 static_cast<vtkIdType>(num_verts),
175 cell_vids.data(),
176 static_cast<vtkIdType>(num_faces),
177 faces_vids.data());
178 break;
179 }
180 default:
181 ugrid->InsertNextCell(vtk_subtype,
182 static_cast<vtkIdType>(num_verts),
183 cell_vids.data());
184 }
185 }//polyhedron
186}
187
188//###################################################################
189/**Uploads vertices and cells to an unstructured grid.*/
191 const std::vector<uint64_t>& vertex_map,
193{
194 const size_t num_verts = cell_face.vertex_ids_.size();
195
196 std::vector<vtkIdType> cell_vids;
197 for (uint64_t vid : cell_face.vertex_ids_)
198 cell_vids.push_back(static_cast<vtkIdType>(vertex_map[vid]));
199
200 if (num_verts == 1)
201 {
202 ugrid->InsertNextCell(VTK_VERTEX,
203 static_cast<vtkIdType>(num_verts),
204 cell_vids.data());
205 }
206 if (num_verts == 2)
207 {
208 ugrid->InsertNextCell(VTK_LINE,
209 static_cast<vtkIdType>(num_verts),
210 cell_vids.data());
211 }
212 if (num_verts >= 3)
213 {
214 int vtk_subtype;
215 switch (num_verts)
216 {
217 case 3: vtk_subtype = VTK_TRIANGLE; break;
218 case 4: vtk_subtype = VTK_QUAD; break;
219 default: vtk_subtype = VTK_POLYGON; break;
220 }
221
222 ugrid->InsertNextCell(vtk_subtype,
223 static_cast<vtkIdType>(num_verts),
224 cell_vids.data());
225 }
226}
std::vector< uint64_t > vertex_ids_
Definition: cell.h:38
std::vector< CellFace > faces_
Definition: cell.h:82
CellType Type() const
Definition: cell.h:98
CellType SubType() const
Definition: cell.h:99
std::vector< uint64_t > vertex_ids_
Definition: cell.h:81
void UploadFaceGeometry(const chi_mesh::CellFace &cell_face, const std::vector< uint64_t > &vertex_map, vtkNew< vtkUnstructuredGrid > &ugrid)
void UploadCellGeometryContinuous(const chi_mesh::Cell &cell, const std::vector< uint64_t > &vertex_map, vtkNew< vtkUnstructuredGrid > &ugrid)
void UploadCellGeometryDiscontinuous(const chi_mesh::MeshContinuum &grid, const chi_mesh::Cell &cell, int64_t &node_counter, vtkNew< vtkPoints > &points, vtkNew< vtkUnstructuredGrid > &ugrid)