Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/ElmerGUI/Application/vtkpost/featureedge.cpp
3203 views
1
/*****************************************************************************
2
* *
3
* Elmer, A Finite Element Software for Multiphysical Problems *
4
* *
5
* Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland *
6
* *
7
* This program is free software; you can redistribute it and/or *
8
* modify it under the terms of the GNU General Public License *
9
* as published by the Free Software Foundation; either version 2 *
10
* of the License, or (at your option) any later version. *
11
* *
12
* This program 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 this program (in file fem/GPL-2); if not, write to the *
19
* Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, *
20
* Boston, MA 02110-1301, USA. *
21
* *
22
*****************************************************************************/
23
24
/*****************************************************************************
25
* *
26
* ElmerGUI featureedge *
27
* *
28
*****************************************************************************
29
* *
30
* Authors: Mikko Lyly, Juha Ruokolainen and Peter R�back *
31
* Email: [email protected] *
32
* Web: http://www.csc.fi/elmer *
33
* Address: CSC - IT Center for Science Ltd. *
34
* Keilaranta 14 *
35
* 02101 Espoo, Finland *
36
* *
37
* Original Date: 15 Mar 2008 *
38
* *
39
*****************************************************************************/
40
41
#include <QtGui>
42
#include <iostream>
43
#include "vtkpost.h"
44
#include "featureedge.h"
45
#include "preferences.h"
46
47
#include <vtkGeometryFilter.h>
48
#include <vtkFeatureEdges.h>
49
#include <vtkPolyDataMapper.h>
50
#include <vtkProperty.h>
51
#include <vtkTubeFilter.h>
52
#include <vtkClipPolyData.h>
53
#include <vtkPlane.h>
54
55
using namespace std;
56
57
FeatureEdge::FeatureEdge(QWidget *parent)
58
: QDialog(parent)
59
{
60
ui.setupUi(this);
61
62
setWindowTitle("Feature edges");
63
setWindowIcon(QIcon(":/icons/Mesh3D.png"));
64
65
actor = vtkActor::New();
66
}
67
68
FeatureEdge::~FeatureEdge()
69
{
70
actor->Delete();
71
}
72
73
/*
74
// The original draw() function which draws all the groups by one FeatureEdge instance
75
// using one vtkUnstructuredGrid. This ends up with boundary of two groups not drawn.
76
void FeatureEdge::draw(VtkPost* vtkPost, Preferences* preferences)
77
{
78
bool useSurfaceGrid = preferences->ui.surfaceRButton->isChecked();
79
int featureAngle = preferences->ui.angleSpin->value();
80
int lineWidth = preferences->ui.lineWidthSpin->value();
81
bool useTubeFilter = preferences->ui.featureEdgeTubes->isChecked();
82
int tubeQuality = preferences->ui.featureEdgeTubeQuality->value();
83
int radius = preferences->ui.featureEdgeTubeRadius->value();
84
bool useClip = preferences->ui.featureEdgesClip->isChecked();
85
bool boundaryEdges = preferences->ui.drawBoundaryEdges->isChecked();
86
useClip |= vtkPost->GetClipAll();
87
88
vtkUnstructuredGrid* grid = NULL;
89
90
if(useSurfaceGrid) {
91
grid = vtkPost->GetSurfaceGrid();
92
} else {
93
grid = vtkPost->GetVolumeGrid();
94
}
95
96
if(!grid) return;
97
98
if(grid->GetNumberOfCells() < 1) return;
99
100
// Convert from vtkUnstructuredGrid to vtkPolyData:
101
vtkGeometryFilter* filter = vtkGeometryFilter::New();
102
#if VTK_MAJOR_VERSION <= 5
103
filter->SetInput(grid);
104
#else
105
filter->SetInputData(grid);
106
#endif
107
// filter->GetOutput()->ReleaseDataFlagOn();
108
109
vtkFeatureEdges* edges = vtkFeatureEdges::New();
110
edges->SetInputConnection(filter->GetOutputPort());
111
edges->SetFeatureAngle(featureAngle);
112
edges->NonManifoldEdgesOn();
113
edges->ManifoldEdgesOn();
114
if(boundaryEdges) {
115
edges->BoundaryEdgesOn();
116
} else {
117
edges->BoundaryEdgesOff();
118
}
119
120
vtkTubeFilter* tubes = vtkTubeFilter::New();
121
if(useTubeFilter) {
122
double r = vtkPost->GetLength() * radius / 2000.0;
123
tubes->SetInputConnection(edges->GetOutputPort());
124
tubes->SetNumberOfSides(tubeQuality);
125
tubes->SetRadius(r);
126
}
127
128
vtkClipPolyData* clipper = vtkClipPolyData::New();
129
if(useClip) {
130
if(useTubeFilter) {
131
clipper->SetInputConnection(tubes->GetOutputPort());
132
} else {
133
clipper->SetInputConnection(edges->GetOutputPort());
134
}
135
clipper->SetClipFunction(vtkPost->GetClipPlane());
136
clipper->GenerateClippedOutputOn();
137
}
138
139
vtkPolyDataMapper* mapper = vtkPolyDataMapper::New();
140
if(useClip) {
141
mapper->SetInputConnection(clipper->GetOutputPort());
142
} else {
143
if(useTubeFilter) {
144
mapper->SetInputConnection(tubes->GetOutputPort());
145
} else {
146
mapper->SetInputConnection(edges->GetOutputPort());
147
}
148
}
149
mapper->ScalarVisibilityOff();
150
mapper->SetResolveCoincidentTopologyToPolygonOffset();
151
// mapper->ImmediateModeRenderingOn();
152
153
vtkPost->GetFeatureEdgeActor()->GetProperty()->SetLineWidth(lineWidth);
154
vtkPost->GetFeatureEdgeActor()->GetProperty()->SetColor(0, 0, 0);
155
vtkPost->GetFeatureEdgeActor()->SetMapper(mapper);
156
157
158
mapper->Delete();
159
clipper->Delete();
160
tubes->Delete();
161
edges->Delete();
162
filter->Delete();
163
}
164
*/
165
166
167
// The new draw() function which draws one group by one FeatureEdge instance
168
// using the specified vtkUnstructuredGrid to draw boundary of two groups.
169
void FeatureEdge::draw(VtkPost* vtkPost, Preferences* preferences, vtkUnstructuredGrid* grid)
170
{
171
bool useSurfaceGrid = preferences->ui.surfaceRButton->isChecked();
172
int featureAngle = preferences->ui.angleSpin->value();
173
int lineWidth = preferences->ui.lineWidthSpin->value();
174
bool useTubeFilter = preferences->ui.featureEdgeTubes->isChecked();
175
int tubeQuality = preferences->ui.featureEdgeTubeQuality->value();
176
int radius = preferences->ui.featureEdgeTubeRadius->value();
177
bool useClip = preferences->ui.featureEdgesClip->isChecked();
178
bool boundaryEdges = preferences->ui.drawBoundaryEdges->isChecked();
179
useClip |= vtkPost->GetClipAll();
180
QColor color = preferences->getFeatureEdgeColor();
181
182
183
if(!grid) return;
184
185
if(grid->GetNumberOfCells() < 1) return;
186
187
// Convert from vtkUnstructuredGrid to vtkPolyData:
188
vtkGeometryFilter* filter = vtkGeometryFilter::New();
189
#if VTK_MAJOR_VERSION <= 5
190
filter->SetInput(grid);
191
#else
192
filter->SetInputData(grid);
193
#endif
194
// filter->GetOutput()->ReleaseDataFlagOn();
195
196
vtkFeatureEdges* edges = vtkFeatureEdges::New();
197
edges->SetInputConnection(filter->GetOutputPort());
198
edges->SetFeatureAngle(featureAngle);
199
edges->NonManifoldEdgesOn();
200
edges->ManifoldEdgesOn();
201
if(boundaryEdges) {
202
edges->BoundaryEdgesOn();
203
} else {
204
edges->BoundaryEdgesOff();
205
}
206
207
vtkTubeFilter* tubes = vtkTubeFilter::New();
208
if(useTubeFilter) {
209
double r = vtkPost->GetLength() * radius / 2000.0;
210
tubes->SetInputConnection(edges->GetOutputPort());
211
tubes->SetNumberOfSides(tubeQuality);
212
tubes->SetRadius(r);
213
}
214
215
vtkClipPolyData* clipper = vtkClipPolyData::New();
216
if(useClip) {
217
if(useTubeFilter) {
218
clipper->SetInputConnection(tubes->GetOutputPort());
219
} else {
220
clipper->SetInputConnection(edges->GetOutputPort());
221
}
222
clipper->SetClipFunction(vtkPost->GetClipPlane());
223
clipper->GenerateClippedOutputOn();
224
}
225
226
vtkPolyDataMapper* mapper = vtkPolyDataMapper::New();
227
if(useClip) {
228
mapper->SetInputConnection(clipper->GetOutputPort());
229
} else {
230
if(useTubeFilter) {
231
mapper->SetInputConnection(tubes->GetOutputPort());
232
} else {
233
mapper->SetInputConnection(edges->GetOutputPort());
234
}
235
}
236
mapper->ScalarVisibilityOff();
237
mapper->SetResolveCoincidentTopologyToPolygonOffset();
238
// mapper->ImmediateModeRenderingOn();
239
240
actor->GetProperty()->SetLineWidth(lineWidth);
241
actor->GetProperty()->SetColor(color.redF(), color.greenF(), color.blueF());
242
actor->SetMapper(mapper);
243
244
mapper->Delete();
245
clipper->Delete();
246
tubes->Delete();
247
edges->Delete();
248
filter->Delete();
249
}
250
251
void FeatureEdge::removeActorFrom(vtkRenderer* renderer){
252
renderer->RemoveActor(actor);
253
}
254
void FeatureEdge::addActorTo(vtkRenderer* renderer){
255
renderer->AddActor(actor);
256
}
257