updated to modern VTK
[engrid-github.git] / src / libengrid / stitchholes.cpp
blob2a0bec20f2498515d8c217171cad9b7da6a3e606
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 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22 #include "stitchholes.h"
23 #include "guimainwindow.h"
24 #include "deletestraynodes.h"
26 #include <vtkCellArray.h>
27 #include <vtkPolyDataWriter.h>
29 StitchHoles::StitchHoles(int bc)
31 m_Bc = bc;
32 m_Cad = GuiMainWindow::pointer()->getCadInterface(bc);
35 QList<vtkIdType> StitchHoles::getNextHole()
37 QList<vtkIdType> loop_nodes;
38 EG_VTKDCC(vtkIntArray, cell_code, m_Grid, "cell_code");
40 // find the first two nodes of the hole
42 vtkIdType id_node1 = -1;
43 vtkIdType id_node2 = -1;
44 EG_FORALL_CELLS (id_cell, m_Grid) {
45 if (isSurface(id_cell, m_Grid) && cell_code->GetValue(id_cell) == m_Bc) {
46 for (int i = 0; i < m_Part.c2cGSize(id_cell); ++i) {
47 if (m_Part.c2cGG(id_cell, i) == -1) {
48 QList<vtkIdType> pts;
49 getPointsOfCell(m_Grid, id_cell, pts);
50 pts << pts.first();
51 id_node1 = pts[i+1];
52 id_node2 = pts[i];
53 break;
56 if (id_node1 != -1) {
57 break;
61 if (id_node1 == -1) {
62 return loop_nodes;
65 // create node loop around hole
67 vtkIdType id_start = id_node1;
68 loop_nodes << id_node1;
69 vec3_t x0;
70 m_Grid->GetPoint(id_node1, x0.data());
71 while (id_node2 != id_start && loop_nodes.size() < m_Grid->GetNumberOfPoints()) {
72 loop_nodes << id_node2;
73 vec3_t x;
74 m_Grid->GetPoint(id_node2, x.data());
75 x0 += x;
76 bool found = false;
77 for (int i = 0; i < m_Part.n2nGSize(id_node2); ++i) {
78 vtkIdType id_node3 = m_Part.n2nGG(id_node2, i);
79 if (id_node3 != id_node1) {
80 QList<vtkIdType> edge_faces;
81 m_Part.getEdgeFaces(id_node2, id_node3, edge_faces);
82 if (edge_faces.size() == 1) {
83 found = true;
84 id_node1 = id_node2;
85 id_node2 = id_node3;
86 break;
90 if (!found) {
91 loop_nodes.clear();
92 return loop_nodes;
95 x0 *= 1.0/loop_nodes.size();
96 m_Cad->snap(x0);
97 vec3_t n = m_Cad->getLastNormal();
98 setOrigin(x0);
99 setNormal(n);
100 setupTransformation();
101 return loop_nodes;
104 vec3_t StitchHoles::transformFromPlane(vec3_t x)
106 vec3_t x_best = x;
107 double dist_min = EG_LARGE_REAL;
108 for (int i = 0; i < m_X2.size(); ++i) {
109 double dist = (x - m_X2[i]).abs();
110 if (dist < dist_min) {
111 dist_min = dist;
112 x_best = m_X3[i];
115 return x_best;
118 void StitchHoles::stitchHole(QList<vtkIdType> loop_nodes)
120 EG_VTKSP(vtkPolyData, edge_pdata);
121 EG_VTKSP(vtkPoints, points);
122 EG_VTKSP(vtkCellArray, polys);
124 m_X2.clear();
125 m_X3.clear();
127 for (int i = 0; i < loop_nodes.size(); ++i) {
128 vec3_t x;
129 m_Grid->GetPoint(loop_nodes[i], x.data());
130 m_X3 << x;
131 x = toPlane(x);
132 m_X2 << x;
133 points->InsertNextPoint(x.data());
135 EG_VTKSP(vtkIdList, pts);
136 pts->SetNumberOfIds(loop_nodes.size());
137 for (vtkIdType i = 0; i < pts->GetNumberOfIds(); ++i) {
138 pts->SetId(i, i);
140 polys->InsertNextCell(pts);
141 edge_pdata->SetPoints(points);
142 edge_pdata->SetPolys(polys);
143 EG_VTKSP(vtkUnstructuredGrid, tri_grid);
144 triangulate(edge_pdata, tri_grid, m_Bc);
145 for (vtkIdType id_node = 0; id_node < tri_grid->GetNumberOfPoints(); ++id_node) {
146 vec3_t x;
147 tri_grid->GetPoint(id_node, x.data());
148 x = transformFromPlane(x);
149 tri_grid->GetPoints()->SetPoint(id_node, x.data());
153 EG_FORALL_CELLS(id_cell, tri_grid) {
154 if (cellNormal(tri_grid, id_cell)*m_N < 0) {
155 EG_GET_CELL(id_cell, tri_grid);
156 QVector<vtkIdType> nodes(num_pts);
157 for (vtkIdType j = 0; j < num_pts; ++j) {
158 nodes[j] = pts[j];
160 for (vtkIdType j = 0; j < num_pts; ++j) {
161 pts[num_pts - j - 1] = nodes[j];
166 MeshPartition tri_part(tri_grid, true);
167 m_Part.addPartition(tri_part);
169 DeleteStrayNodes del_stray;
170 del_stray.setGrid(m_Grid);
171 del_stray.setAllCells();
172 del_stray();
174 ++m_Counter;
175 QString counter_txt;
176 counter_txt.setNum(m_Counter);
177 counter_txt = counter_txt.rightJustified(3, '0');
180 void StitchHoles::operate()
182 m_Counter = 0;
183 bool hole_found = false;
184 cout << "stitching holes of \"" << qPrintable(GuiMainWindow::pointer()->getBC(m_Bc).getName()) << "\"" << endl;
185 int count = 0;
186 m_Part.trackGrid(m_Grid);
187 do {
188 QList<vtkIdType> loop_nodes = getNextHole();
189 hole_found = loop_nodes.size() > 0;
190 if (hole_found) {
191 stitchHole(loop_nodes);
192 ++count;
193 cout << " " << count << ". hole stiched" << endl;
195 } while (hole_found && count < 50);
196 cout << " finished" << endl;