1 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + This file is part of enGrid. +
5 // + Copyright 2008-2014 enGits GmbH +
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. +
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. +
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/>. +
20 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 #include "stlreader.h"
22 #include "correctsurfaceorientation.h"
24 #include "vtkEgPolyDataToUnstructuredGridFilter.h"
25 #include <vtkSTLReader.h>
26 #include <vtkCleanPolyData.h>
27 #include <vtkFeatureEdges.h>
30 #include <QInputDialog>
32 #include "guimainwindow.h"
33 #include "fixcadgeometry.h"
35 StlReader::StlReader()
37 setFormat("STL files(*.stl *.STL)");
39 m_FileNameSet
= false;
40 m_MaxNumCleanIter
= 20;
43 void StlReader::setFileName(QString file_name
)
45 m_FileName
= file_name
;
46 cout
<< qPrintable(m_FileName
) << endl
;
50 void StlReader::operate()
52 QFileInfo
file_info(GuiMainWindow::pointer()->getFilename());
55 file_name
= m_FileName
;
57 readInputFileName(file_info
.completeBaseName() + ".stl");
59 file_name
= getFileName();
65 EG_VTKSP(vtkSTLReader
, stl
);
67 stl
->SetFileName(qPrintable(file_name
));
69 EG_VTKSP(vtkPolyData
, poly
);
70 poly
->DeepCopy(stl
->GetOutput());
73 for (vtkIdType cellId
= 0; cellId
< poly
->GetNumberOfCells(); ++cellId
) {
75 poly
->GetCellPoints(cellId
, Npts
, pts
);
76 for (int i
= 0; i
< Npts
; ++i
) {
78 poly
->GetPoints()->GetPoint(pts
[i
], x1
.data());
80 poly
->GetPoints()->GetPoint(pts
[0], x2
.data());
82 poly
->GetPoints()->GetPoint(pts
[i
+1], x2
.data());
84 L
= min(L
, (x1
-x2
).abs());
87 if (m_Tolerance
< 0) {
88 m_Tolerance
= QInputDialog::getText(NULL
, "enter STL tolerance", "tolerance", QLineEdit::Normal
, "1e-10").toDouble();
90 cout
<< "cleaning STL geometry:" << endl
;
91 EG_VTKSP(vtkCleanPolyData
, poly_clean
);
92 EG_VTKSP(vtkFeatureEdges
, topo_check
);
94 poly
->GetBounds(bounds
);
95 poly_clean
->ToleranceIsAbsoluteOn();
96 poly_clean
->ConvertLinesToPointsOn();
97 poly_clean
->ConvertPolysToLinesOn();
98 poly_clean
->SetInputData(poly
);
99 topo_check
->SetInputConnection(poly_clean
->GetOutputPort());
100 topo_check
->BoundaryEdgesOn();
101 topo_check
->ManifoldEdgesOff();
102 topo_check
->FeatureEdgesOff();
103 topo_check
->NonManifoldEdgesOn();
108 cout
<< " tolerance = " << m_Tolerance
<< endl
;
109 poly_clean
->SetAbsoluteTolerance(m_Tolerance
);
110 topo_check
->Update();
112 check_passed
= topo_check
->GetOutput()->GetNumberOfPoints() == 0;
113 } while (m_Tolerance
< 1 && !check_passed
&& count
< m_MaxNumCleanIter
);
115 cout
<< "The STL geometry seems to be clean." << endl
;
117 cout
<< "The STL geometry could not be cleaned." << endl
;
119 // with a tolerance of " << 0.5*L << endl;
120 EG_VTKSP(vtkEgPolyDataToUnstructuredGridFilter
, poly2ugrid
);
121 poly2ugrid
->SetInputConnection(poly_clean
->GetOutputPort());
122 poly2ugrid
->Update();
124 allocateGrid(m_Grid
, poly2ugrid
->GetOutput()->GetNumberOfCells(), poly2ugrid
->GetOutput()->GetNumberOfPoints());
125 for (vtkIdType id_node
= 0; id_node
< poly2ugrid
->GetOutput()->GetNumberOfPoints(); ++id_node
) {
127 poly2ugrid
->GetOutput()->GetPoints()->GetPoint(id_node
, x
.data());
128 m_Grid
->GetPoints()->SetPoint(id_node
, x
.data());
130 for (vtkIdType id_cell
= 0; id_cell
< poly2ugrid
->GetOutput()->GetNumberOfCells(); ++id_cell
) {
131 vtkIdType N_pts
, *pts
;
132 vtkIdType type_cell
= poly2ugrid
->GetOutput()->GetCellType(id_cell
);
133 poly2ugrid
->GetOutput()->GetCellPoints(id_cell
, N_pts
, pts
);
134 m_Grid
->InsertNextCell(type_cell
, N_pts
, pts
);
137 EG_VTKDCC(vtkIntArray
, bc
, m_Grid
, "cell_code");
138 EG_VTKDCC(vtkIntArray
, orgdir
, m_Grid
, "cell_orgdir");
139 EG_VTKDCC(vtkIntArray
, voldir
, m_Grid
, "cell_voldir");
140 EG_VTKDCC(vtkIntArray
, curdir
, m_Grid
, "cell_curdir");
141 for (vtkIdType id_cell
= 0; id_cell
< m_Grid
->GetNumberOfCells(); ++id_cell
) {
142 bc
->SetValue(id_cell
, 1);
143 orgdir
->SetValue(id_cell
, 0);
144 voldir
->SetValue(id_cell
, 0);
145 curdir
->SetValue(id_cell
, 0);
148 CorrectSurfaceOrientation corr_surf
;
149 corr_surf
.setGrid(m_Grid
);
151 FixCadGeometry cad_fix
;
152 cad_fix
.setGrid(m_Grid
);