updated to modern VTK
[engrid-github.git] / src / libengrid / vtkEgExtractVolumeCells.cxx
blob1539fb51a8339cb55292085db935e80b5b2dde7d
1 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 // + +
3 // + This file is part of enGrid. +
4 // + +
5 // + Copyright 2008-2014 enGits GmbH +
6 // + +
7 // + enGrid is free software: you can redistribute it and/or modify +
8 // + it under the terms of the GNU General Public License as published by +
9 // + the Free Software Foundation, either version 3 of the License, or +
10 // + (at your option) any later version. +
11 // + +
12 // + enGrid is distributed in the hope that it will be useful, +
13 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
14 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
15 // + GNU General Public License for more details. +
16 // + +
17 // + You should have received a copy of the GNU General Public License +
18 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
19 // + +
20 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 #include "vtkEgExtractVolumeCells.h"
23 #include <vtkIdList.h>
24 #include <vtkSmartPointer.h>
25 #include <vtkType.h>
27 vtkStandardNewMacro(vtkEgExtractVolumeCells)
29 vtkEgExtractVolumeCells::vtkEgExtractVolumeCells()
31 //SetClippingOff uses an IF to check the value first, which may have unexpected behavior in some OSes.
32 m_Clip = false; SetClippingOff();
33 SetX(vec3_t(0,0,0));
34 SetN(vec3_t(1,0,0));
35 m_ExtrTetras = true;
36 m_ExtrPyramids = true;
37 m_ExtrWedges = true;
38 m_ExtrHexes = true;
39 m_ExtrPolys = true;
42 void vtkEgExtractVolumeCells::SetX(vec3_t x)
44 m_X = x;
45 Modified();
48 void vtkEgExtractVolumeCells::SetN(vec3_t n)
50 m_N = n;
51 Modified();
54 void vtkEgExtractVolumeCells::SetClippingOn()
56 if (!m_Clip) {
57 m_Clip = true;
58 Modified();
62 void vtkEgExtractVolumeCells::SetClippingOff()
64 if (m_Clip) {
65 m_Clip = false;
66 Modified();
71 void vtkEgExtractVolumeCells::SetAllOn()
73 SetTetrasOn();
74 SetPyramidsOn();
75 SetWedgesOn();
76 SetHexesOn();
77 SetPolysOn();
80 void vtkEgExtractVolumeCells::SetAllOff()
82 SetTetrasOff();
83 SetPyramidsOff();
84 SetWedgesOff();
85 SetHexesOff();
86 SetPolysOff();
89 void vtkEgExtractVolumeCells::SetTetrasOn()
91 if (!m_ExtrTetras) {
92 m_ExtrTetras = true;
93 Modified();
97 void vtkEgExtractVolumeCells::SetTetrasOff()
99 if (m_ExtrTetras) {
100 m_ExtrTetras = false;
101 Modified();
105 void vtkEgExtractVolumeCells::SetPyramidsOn()
107 if (!m_ExtrPyramids) {
108 m_ExtrPyramids = true;
109 Modified();
113 void vtkEgExtractVolumeCells::SetPyramidsOff()
115 if (m_ExtrPyramids) {
116 m_ExtrPyramids = false;
117 Modified();
121 void vtkEgExtractVolumeCells::SetWedgesOn()
123 if (!m_ExtrWedges) {
124 m_ExtrWedges = true;
125 Modified();
129 void vtkEgExtractVolumeCells::SetWedgesOff()
131 if (m_ExtrWedges) {
132 m_ExtrWedges = false;
133 Modified();
137 void vtkEgExtractVolumeCells::SetHexesOn()
139 if (!m_ExtrHexes) {
140 m_ExtrHexes = true;
141 Modified();
145 void vtkEgExtractVolumeCells::SetHexesOff()
147 if (m_ExtrHexes) {
148 m_ExtrHexes = false;
149 Modified();
153 void vtkEgExtractVolumeCells::SetPolysOn()
155 if (!m_ExtrPolys) {
156 m_ExtrPolys = true;
157 Modified();
161 void vtkEgExtractVolumeCells::SetPolysOff()
163 if (m_ExtrPolys) {
164 m_ExtrPolys = false;
165 Modified();
169 void vtkEgExtractVolumeCells::Setx(double x)
171 m_X[0] = x;
172 Modified();
175 void vtkEgExtractVolumeCells::Sety(double y)
177 m_X[1] = y;
178 Modified();
181 void vtkEgExtractVolumeCells::Setz(double z)
183 m_X[2] = z;
184 Modified();
187 void vtkEgExtractVolumeCells::Setnx(double nx)
189 m_N[0] = nx;
190 Modified();
193 void vtkEgExtractVolumeCells::Setny(double ny)
195 m_N[1] = ny;
196 Modified();
199 void vtkEgExtractVolumeCells::Setnz(double nz)
201 m_N[2] = nz;
202 Modified();
206 void vtkEgExtractVolumeCells::ExecuteEg()
208 QSet<vtkIdType> ex_cells;
209 for (vtkIdType id_cell = 0; id_cell < m_Input->GetNumberOfCells(); ++id_cell) {
210 if (isVolume(id_cell, m_Input)) {
211 bool select = true;
212 vtkIdType type_cell = m_Input->GetCellType(id_cell);
213 if (!m_ExtrTetras && type_cell == VTK_TETRA) {
214 select = false;
216 if (!m_ExtrPyramids && type_cell == VTK_PYRAMID) {
217 select = false;
219 if (!m_ExtrWedges && type_cell == VTK_WEDGE) {
220 select = false;
222 if (!m_ExtrHexes && type_cell == VTK_HEXAHEDRON) {
223 select = false;
225 if (!m_ExtrPolys && type_cell == VTK_POLYHEDRON) {
226 select = false;
228 if (m_Clip && select) {
229 QList<vtkIdType> pts;
230 getPointsOfCell(m_Input, id_cell, pts);
231 foreach (vtkIdType id_node, pts) {
232 vec3_t x;
233 m_Input->GetPoints()->GetPoint(id_node, x.data());
234 if ((x - m_X)*m_N < 0) {
235 select = false;
236 break;
240 if (select) {
241 ex_cells.insert(id_cell);
245 QVector<vtkIdType> cells(ex_cells.size());
246 qCopy(ex_cells.begin(), ex_cells.end(), cells.begin());
247 QVector<vtkIdType> nodes;
248 QVector<int> _nodes;
249 getNodesFromCells(cells, nodes, m_Input);
250 createNodeMapping(nodes, _nodes, m_Input);
251 allocateGrid(m_Output, cells.size(), nodes.size());
252 foreach(vtkIdType id_node, nodes) {
253 vec3_t x;
254 m_Input->GetPoints()->GetPoint(id_node, x.data());
255 m_Output->GetPoints()->SetPoint(_nodes[id_node], x.data());
257 foreach(vtkIdType id_cell, cells) {
258 vtkSmartPointer<vtkIdList> stream = vtkSmartPointer<vtkIdList>::New();
259 vtkSmartPointer<vtkIdList> new_stream = vtkSmartPointer<vtkIdList>::New();
260 vtkIdType type_cell = m_Input->GetCellType(id_cell);
261 m_Input->GetFaceStream(id_cell, stream);
262 new_stream->SetNumberOfIds(stream->GetNumberOfIds());
263 if (type_cell == VTK_POLYHEDRON) {
264 vtkIdType num = stream->GetId(0);
266 new_stream->SetId(0, stream->GetId(0));
267 vtkIdType id = 1;
268 for (int i = 0; i < num; ++i) {
269 vtkIdType num_pts = stream->GetId(id);
270 new_stream->SetId(id, stream->GetId(id));
271 ++id;
272 for (int j = 0; j < num_pts; ++j) {
273 new_stream->SetId(id, _nodes[stream->GetId(id)]);
274 ++id;
278 } else {
279 for (vtkIdType i = 0; i < stream->GetNumberOfIds(); ++i) {
280 new_stream->SetId(i, _nodes[stream->GetId(i)]);
284 vtkIdType id_new_cell = m_Output->InsertNextCell(type_cell, new_stream);
285 copyCellData(m_Input, id_cell, m_Output, id_new_cell);