Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/ElmerGUI/Application/vtkpost/surface.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 surface *
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 <QColorDialog>
43
#include <iostream>
44
#include "epmesh.h"
45
#include "vtkpost.h"
46
#include "surface.h"
47
#include "timestep.h"
48
49
#include <vtkUnstructuredGrid.h>
50
#include <vtkPointData.h>
51
#include <vtkFloatArray.h>
52
#include <vtkGeometryFilter.h>
53
#include <vtkClipPolyData.h>
54
#include <vtkPlane.h>
55
#include <vtkPolyDataNormals.h>
56
#include <vtkDataSetMapper.h>
57
#include <vtkProperty.h>
58
#include <vtkActor.h>
59
60
using namespace std;
61
62
Surface::Surface(QWidget *parent)
63
: QDialog(parent)
64
{
65
ui.setupUi(this);
66
67
connect(ui.cancelButton, SIGNAL(clicked()), this, SLOT(cancelButtonClicked()));
68
connect(ui.applyButton, SIGNAL(clicked()), this, SLOT(applyButtonClicked()));
69
connect(ui.okButton, SIGNAL(clicked()), this, SLOT(okButtonClicked()));
70
connect(ui.surfaceCombo, SIGNAL(currentIndexChanged(int)), this, SLOT(surfaceSelectionChanged(int)));
71
connect(ui.keepLimits, SIGNAL(stateChanged(int)), this, SLOT(keepLimitsSlot(int)));
72
73
setWindowIcon(QIcon(":/icons/Mesh3D.png"));
74
75
ui.cancelButton->setIcon(QIcon::fromTheme("dialog-error-round"));
76
ui.applyButton->setIcon(QIcon::fromTheme("view-refresh"));
77
ui.okButton->setIcon(QIcon::fromTheme("dialog-accept"));
78
79
setNullColor(Qt::blue);
80
connect(ui.nullColorButton, SIGNAL(clicked()), this, SLOT(nullColorButtonClicked()));
81
}
82
83
Surface::~Surface()
84
{
85
}
86
87
void Surface::cancelButtonClicked()
88
{
89
emit(hideSurfaceSignal());
90
close();
91
}
92
93
void Surface::applyButtonClicked()
94
{
95
emit(drawSurfaceSignal());
96
}
97
98
void Surface::okButtonClicked()
99
{
100
applyButtonClicked();
101
close();
102
}
103
104
void Surface::populateWidgets(VtkPost* vtkPost)
105
{
106
this->scalarField = vtkPost->GetScalarField();
107
this->scalarFields = vtkPost->GetScalarFields();
108
109
QString name = ui.surfaceCombo->currentText();
110
111
ui.surfaceCombo->clear();
112
113
for(int i = 0; i < scalarFields; i++) {
114
ScalarField *sf = &scalarField[i];
115
QString name = sf->name;
116
ui.surfaceCombo->addItem(sf->name);
117
}
118
119
this->SetFieldName(name);
120
121
surfaceSelectionChanged(ui.surfaceCombo->currentIndex());
122
}
123
124
void Surface::surfaceSelectionChanged(int newIndex)
125
{
126
ScalarField *sf = &this->scalarField[newIndex];
127
if(!ui.keepLimits->isChecked()) {
128
ui.minEdit->setText(QString::number(sf->minVal));
129
ui.maxEdit->setText(QString::number(sf->maxVal));
130
}
131
132
if(ui.surfaceCombo->currentIndex() == 0 ){ // i.e. Null field
133
ui.nullColorLabel->show();
134
ui.nullColorButton->show();
135
ui.minEdit->setEnabled(false);
136
ui.maxEdit->setEnabled(false);
137
ui.minLabel->setEnabled(false);
138
ui.maxLabel->setEnabled(false);
139
ui.keepLimits->setEnabled(false);
140
}else{
141
ui.nullColorLabel->hide();
142
ui.nullColorButton->hide();
143
ui.minEdit->setEnabled(true);
144
ui.maxEdit->setEnabled(true);
145
ui.minLabel->setEnabled(true);
146
ui.maxLabel->setEnabled(true);
147
ui.keepLimits->setEnabled(true);
148
}
149
}
150
151
void Surface::keepLimitsSlot(int state)
152
{
153
if(state == 0)
154
surfaceSelectionChanged(ui.surfaceCombo->currentIndex());
155
}
156
157
void Surface::draw(VtkPost* vtkPost, TimeStep* timeStep)
158
{
159
int surfaceIndex = ui.surfaceCombo->currentIndex();
160
QString surfaceName = ui.surfaceCombo->currentText();
161
double minVal = ui.minEdit->text().toDouble();
162
double maxVal = ui.maxEdit->text().toDouble();
163
bool useNormals = ui.useNormals->isChecked();
164
int featureAngle = ui.featureAngle->value();
165
double opacity = ui.opacitySpin->value() / 100.0;
166
bool useClip = ui.clipPlane->isChecked();
167
useClip |= vtkPost->GetClipAll();
168
169
ScalarField* sf = &scalarField[surfaceIndex];
170
int maxDataSteps = sf->values / vtkPost->NofNodes();
171
int step = timeStep->ui.timeStep->value();
172
if(step > maxDataSteps) step = maxDataSteps;
173
if(step > timeStep->maxSteps) step = timeStep->maxSteps;
174
int offset = vtkPost->NofNodes() * (step-1);
175
176
// Scalars:
177
//---------
178
vtkPost->GetSurfaceGrid()->GetPointData()->RemoveArray("Surface");
179
vtkFloatArray* scalars = vtkFloatArray::New();
180
scalars->SetNumberOfComponents(1);
181
scalars->SetNumberOfTuples(vtkPost->NofNodes());
182
scalars->SetName("Surface");
183
for(int i = 0; i < vtkPost->NofNodes(); i++)
184
scalars->SetComponent(i, 0, sf->value[i + offset]);
185
vtkPost->GetSurfaceGrid()->GetPointData()->AddArray(scalars);
186
187
// Convert from vtkUnstructuredGrid to vtkPolyData:
188
//-------------------------------------------------
189
vtkGeometryFilter* filter = vtkGeometryFilter::New();
190
191
#if VTK_MAJOR_VERSION <= 5
192
filter->SetInput(vtkPost->GetSurfaceGrid());
193
#else
194
filter->SetInputData(vtkPost->GetSurfaceGrid());
195
#endif
196
#if VTK_MAJOR_VERSION <= 5
197
filter->GetOutput()->ReleaseDataFlagOn();
198
#else
199
filter->ReleaseDataFlagOn();
200
#endif
201
202
// Apply the clip plane:
203
//-----------------------
204
vtkClipPolyData *clipper = vtkClipPolyData::New();
205
206
if(useClip) {
207
clipper->SetInputConnection(filter->GetOutputPort());
208
clipper->SetClipFunction(vtkPost->GetClipPlane());
209
clipper->GenerateClipScalarsOn();
210
clipper->GenerateClippedOutputOn();
211
}
212
213
// Normals:
214
//---------
215
vtkPolyDataNormals *normals = vtkPolyDataNormals::New();
216
217
if(useNormals) {
218
if(useClip) {
219
normals->SetInputConnection(clipper->GetOutputPort());
220
} else {
221
normals->SetInputConnection(filter->GetOutputPort());
222
}
223
normals->SetFeatureAngle(featureAngle);
224
}
225
226
// Mapper:
227
//--------
228
vtkDataSetMapper *mapper = vtkDataSetMapper::New();
229
230
if(useNormals) {
231
mapper->SetInputConnection(normals->GetOutputPort());
232
} else {
233
if(useClip) {
234
mapper->SetInputConnection(clipper->GetOutputPort());
235
} else {
236
#if VTK_MAJOR_VERSION <= 5
237
mapper->SetInput(vtkPost->GetSurfaceGrid());
238
#else
239
mapper->SetInputData(vtkPost->GetSurfaceGrid());
240
#endif
241
}
242
}
243
244
mapper->SetScalarModeToUsePointFieldData();
245
mapper->SelectColorArray("Surface");
246
mapper->ScalarVisibilityOn();
247
mapper->SetScalarRange(minVal, maxVal);
248
mapper->SetResolveCoincidentTopologyToPolygonOffset();
249
mapper->InterpolateScalarsBeforeMappingOn();
250
//mapper->SetLookupTable(vtkPost->GetCurrentLut());
251
mapper->SetLookupTable(vtkPost->GetLut("Surface"));
252
// mapper->ImmediateModeRenderingOn();
253
if(ui.surfaceCombo->currentIndex() == 0 ){ // i.e. Null field
254
mapper->SetScalarRange(0, 1);
255
double h = nullColor.hueF();
256
double s = nullColor.saturationF();
257
double v = nullColor.valueF();
258
int nColor =128;
259
vtkLookupTable* nullLut = vtkLookupTable::New();
260
nullLut->SetHueRange(h, h);
261
nullLut->SetSaturationRange(s, s);
262
nullLut->SetValueRange(v, v);
263
nullLut->SetNumberOfColors(nColor);
264
nullLut->Build();
265
mapper->SetLookupTable(nullLut);
266
nullLut->Delete();
267
}
268
269
// Actor:
270
//--------
271
vtkPost->GetSurfaceActor()->SetMapper(mapper);
272
vtkPost->GetSurfaceActor()->GetProperty()->SetOpacity(opacity);
273
vtkPost->SetCurrentSurfaceName(sf->name);
274
275
// Clean up:
276
//-----------
277
mapper->Delete();
278
normals->Delete();
279
clipper->Delete();
280
filter->Delete();
281
scalars->Delete();
282
}
283
284
// Public slots:
285
//---------------
286
QString Surface::GetFieldName()
287
{
288
return ui.surfaceCombo->currentText();
289
}
290
291
bool Surface::SetFieldName(QString name)
292
{
293
for(int i = 0; i < ui.surfaceCombo->count(); i++) {
294
if(ui.surfaceCombo->itemText(i) == name) {
295
ui.surfaceCombo->setCurrentIndex(i);
296
return true;
297
}
298
}
299
return false;
300
}
301
302
void Surface::SetMinVal(double f)
303
{
304
ui.minEdit->setText(QString::number(f));
305
}
306
307
void Surface::SetMaxVal(double f)
308
{
309
ui.maxEdit->setText(QString::number(f));
310
}
311
312
void Surface::KeepLimits(bool b)
313
{
314
ui.keepLimits->setChecked(b);
315
}
316
317
void Surface::SetComputeNormals(bool b)
318
{
319
ui.useNormals->setChecked(b);
320
}
321
322
void Surface::SetFeatureAngle(int n)
323
{
324
ui.featureAngle->setValue(n);
325
}
326
327
void Surface::SetOpacity(int n)
328
{
329
ui.opacitySpin->setValue(n);
330
}
331
332
void Surface::SetClipPlane(bool b)
333
{
334
ui.clipPlane->setChecked(b);
335
}
336
337
void Surface::nullColorButtonClicked()
338
{
339
setNullColor(QColorDialog::getColor(nullColor));
340
}
341
342
void Surface::setNullColor(QColor color){
343
if(!color.isValid()) return;
344
345
nullColor = color;
346
347
QPalette plt(ui.nullColorLabel->palette());
348
plt.setColor(QPalette::WindowText, nullColor);
349
ui.nullColorLabel->setPalette(plt);
350
}
351