Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/ElmerGUI/Application/vtkpost/streamline.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 streamline *
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 "streamline.h"
47
#include "timestep.h"
48
49
#include <vtkUnstructuredGrid.h>
50
#include <vtkFloatArray.h>
51
#include <vtkPointData.h>
52
#include <vtkPointSource.h>
53
#include <vtkLineSource.h>
54
55
#if VTK_MAJOR_VERSION <= 7
56
#include <vtkStreamLine.h>
57
#else
58
#include <vtkStreamTracer.h>
59
#endif
60
61
#include <vtkRungeKutta4.h>
62
#include <vtkRibbonFilter.h>
63
#include <vtkPolyDataMapper.h>
64
#include <vtkProperty.h>
65
#include <vtkLookupTable.h>
66
#include <vtkActor.h>
67
68
using namespace std;
69
70
StreamLine::StreamLine(QWidget *parent)
71
: QDialog(parent)
72
{
73
ui.setupUi(this);
74
75
connect(ui.cancelButton, SIGNAL(clicked()), this, SLOT(cancelButtonClicked()));
76
connect(ui.applyButton, SIGNAL(clicked()), this, SLOT(applyButtonClicked()));
77
connect(ui.okButton, SIGNAL(clicked()), this, SLOT(okButtonClicked()));
78
connect(ui.colorCombo, SIGNAL(currentIndexChanged(int)), this, SLOT(colorSelectionChanged(int)));
79
connect(ui.keepLimits, SIGNAL(stateChanged(int)), this, SLOT(keepLimitsSlot(int)));
80
81
setWindowIcon(QIcon(":/icons/Mesh3D.png"));
82
83
ui.cancelButton->setIcon(QIcon::fromTheme("dialog-error-round"));
84
ui.applyButton->setIcon(QIcon::fromTheme("view-refresh"));
85
ui.okButton->setIcon(QIcon::fromTheme("dialog-accept"));
86
87
#if VTK_MAJOR_VERSION >= 8
88
ui.propagLabel->setText(QString("Max propagation:"));
89
ui.stepLabel->setText(QString("Max integ. step:"));
90
ui.integStepLabel->setText(QString("Initial integ. step:"));
91
ui.threadsLabel->setText(QString("Max step:"));
92
ui.threads->setValue(2000);
93
#endif
94
95
setNullColor(Qt::blue);
96
connect(ui.nullColorButton, SIGNAL(clicked()), this, SLOT(nullColorButtonClicked()));
97
}
98
99
StreamLine::~StreamLine()
100
{
101
}
102
103
void StreamLine::applyButtonClicked()
104
{
105
emit(drawStreamLineSignal());
106
}
107
108
void StreamLine::cancelButtonClicked()
109
{
110
emit(hideStreamLineSignal());
111
close();
112
}
113
114
void StreamLine::okButtonClicked()
115
{
116
applyButtonClicked();
117
close();
118
}
119
120
void StreamLine::populateWidgets(VtkPost* vtkPost)
121
{
122
this->scalarField = vtkPost->GetScalarField();
123
this->scalarFields = vtkPost->GetScalarFields();
124
125
ui.vectorCombo->clear();
126
127
int index = -1;
128
for(int i = 0; i < scalarFields; i++) {
129
ScalarField *sf = &scalarField[i];
130
QString name = sf->name;
131
if((index = name.indexOf("_x")) >= 0) {
132
ui.vectorCombo->addItem(name.mid(0, index));
133
}
134
}
135
136
QString name = ui.colorCombo->currentText();
137
138
ui.colorCombo->clear();
139
140
for(int i = 0; i < scalarFields; i++) {
141
ScalarField *sf = &scalarField[i];
142
ui.colorCombo->addItem(sf->name);
143
}
144
145
this->SetColorName(name);
146
147
colorSelectionChanged(ui.colorCombo->currentIndex());
148
}
149
150
void StreamLine::colorSelectionChanged(int newIndex)
151
{
152
ScalarField *sf = &this->scalarField[newIndex];
153
if(!ui.keepLimits->isChecked()) {
154
ui.minVal->setText(QString::number(sf->minVal));
155
ui.maxVal->setText(QString::number(sf->maxVal));
156
}
157
if(ui.colorCombo->currentIndex() == 0 ){ // i.e. Null field
158
ui.nullColorLabel->show();
159
ui.nullColorButton->show();
160
ui.minVal->setEnabled(false);
161
ui.maxVal->setEnabled(false);
162
ui.minLabel->setEnabled(false);
163
ui.maxLabel->setEnabled(false);
164
ui.keepLimits->setEnabled(false);
165
}else{
166
ui.nullColorLabel->hide();
167
ui.nullColorButton->hide();
168
ui.minVal->setEnabled(true);
169
ui.maxVal->setEnabled(true);
170
ui.minLabel->setEnabled(true);
171
ui.maxLabel->setEnabled(true);
172
ui.keepLimits->setEnabled(true);
173
}
174
}
175
176
void StreamLine::keepLimitsSlot(int state)
177
{
178
if(state == 0)
179
colorSelectionChanged(ui.colorCombo->currentIndex());
180
}
181
182
void StreamLine::draw(VtkPost* vtkPost, TimeStep* timeStep)
183
{
184
185
QString vectorName = ui.vectorCombo->currentText();
186
187
if(vectorName.isEmpty()) return;
188
189
int i, j, index = -1;
190
for(i = 0; i < scalarFields; i++) {
191
ScalarField* sf = &scalarField[i];
192
QString name = sf->name;
193
if((j = name.indexOf("_x")) >= 0) {
194
if(vectorName == name.mid(0, j)) {
195
index = i;
196
break;
197
}
198
}
199
}
200
201
if(index < 0) return;
202
203
// Controls:
204
double propagationTime = ui.propagationTime->text().toDouble();
205
double stepLength = ui.stepLength->text().toDouble();
206
double integStepLength = ui.integStepLength->text().toDouble();
207
int threads = ui.threads->value();
208
bool useSurfaceGrid = ui.useSurfaceGrid->isChecked();
209
bool lineSource = ui.lineSource->isChecked();
210
bool sphereSource = ui.sphereSource->isChecked();
211
bool pickSource = ui.pickSource->isChecked();
212
bool forward = ui.forward->isChecked();
213
bool backward = ui.backward->isChecked();
214
215
if(!(forward || backward)) return;
216
217
// Color:
218
int colorIndex = ui.colorCombo->currentIndex();
219
QString colorName = ui.colorCombo->currentText();
220
double minVal = ui.minVal->text().toDouble();
221
double maxVal = ui.maxVal->text().toDouble();
222
223
// Appearance:
224
bool drawRibbon = ui.drawRibbon->isChecked();
225
int ribbonWidth = ui.ribbonWidth->value();
226
int lineWidth = ui.lineWidth->text().toInt();
227
228
// Line source:
229
double startX = ui.startX->text().toDouble();
230
double startY = ui.startY->text().toDouble();
231
double startZ = ui.startZ->text().toDouble();
232
double endX = ui.endX->text().toDouble();
233
double endY = ui.endY->text().toDouble();
234
double endZ = ui.endZ->text().toDouble();
235
int lines = ui.lines->value();
236
bool drawRake = ui.rake->isChecked();
237
int rakeWidth = ui.rakeWidth->value();
238
239
ScalarField* sf_x = &scalarField[index + 0];
240
ScalarField* sf_y = &scalarField[index + 1];
241
ScalarField* sf_z = &scalarField[index + 2];
242
int maxDataStepVector = sf_x->values / vtkPost->NofNodes();
243
int step = timeStep->ui.timeStep->value();
244
if(step > maxDataStepVector) step = maxDataStepVector;
245
if(step > timeStep->maxSteps) step = timeStep->maxSteps;
246
int vectorOffset = vtkPost->NofNodes() * (step-1);
247
248
ScalarField* sf = &scalarField[colorIndex];
249
int maxDataStepColor = sf->values / vtkPost->NofNodes();
250
step = timeStep->ui.timeStep->value();
251
if(step > maxDataStepColor) step = maxDataStepColor;
252
if(step > timeStep->maxSteps) step = timeStep->maxSteps;
253
int colorOffset = vtkPost->NofNodes() * (step-1);
254
255
// Sphere source:
256
double centerX = ui.centerX->text().toDouble();
257
double centerY = ui.centerY->text().toDouble();
258
double centerZ = ui.centerZ->text().toDouble();
259
double radius = ui.radius->text().toDouble();
260
int points = ui.points->value();
261
262
// Pick source:
263
double* currentPickPosition = vtkPost->GetCurrentPickPosition();
264
double pickX = currentPickPosition[0];
265
double pickY = currentPickPosition[1];
266
double pickZ = currentPickPosition[2];
267
268
// Choose the grid:
269
//------------------
270
vtkUnstructuredGrid* grid = NULL;
271
if(useSurfaceGrid)
272
grid = vtkPost->GetSurfaceGrid();
273
else
274
grid = vtkPost->GetVolumeGrid();
275
276
if(!grid) return;
277
278
if(grid->GetNumberOfCells() < 1) return;
279
280
// Vector data:
281
//-------------
282
grid->GetPointData()->RemoveArray("VectorData");
283
vtkFloatArray* vectorData = vtkFloatArray::New();
284
vectorData->SetNumberOfComponents(3);
285
vectorData->SetNumberOfTuples(vtkPost->NofNodes());
286
vectorData->SetName("VectorData");
287
for(int i = 0; i < vtkPost->NofNodes(); i++) {
288
double val_x = sf_x->value[i + vectorOffset];
289
double val_y = sf_y->value[i + vectorOffset];
290
double val_z = sf_z->value[i + vectorOffset];
291
vectorData->SetComponent(i,0,val_x);
292
vectorData->SetComponent(i,1,val_y);
293
vectorData->SetComponent(i,2,val_z);
294
}
295
grid->GetPointData()->AddArray(vectorData);
296
297
// Color data:
298
//-------------
299
grid->GetPointData()->RemoveArray("StreamLineColor");
300
vtkFloatArray* vectorColor = vtkFloatArray::New();
301
vectorColor->SetNumberOfComponents(1);
302
vectorColor->SetNumberOfTuples(vtkPost->NofNodes());
303
vectorColor->SetName("StreamLineColor");
304
for(int i = 0; i < vtkPost->NofNodes(); i++)
305
vectorColor->SetComponent(i, 0, sf->value[i + colorOffset]);
306
grid->GetPointData()->AddArray(vectorColor);
307
308
// Generate stream lines:
309
//-----------------------
310
grid->GetPointData()->SetActiveVectors("VectorData");
311
grid->GetPointData()->SetActiveScalars("StreamLineColor");
312
vtkPointSource* point = vtkPointSource::New();
313
vtkLineSource* line = vtkLineSource::New();
314
if(lineSource) {
315
line->SetPoint1(startX, startY, startZ);
316
line->SetPoint2(endX, endY, endZ);
317
line->SetResolution(lines);
318
} else {
319
if(sphereSource) {
320
point->SetCenter(centerX, centerY, centerZ);
321
point->SetRadius(radius);
322
point->SetNumberOfPoints(points);
323
point->SetDistributionToUniform();
324
} else {
325
point->SetCenter(pickX, pickY, pickZ);
326
point->SetRadius(0.0);
327
point->SetNumberOfPoints(1);
328
}
329
}
330
331
#if VTK_MAJOR_VERSION <= 7
332
vtkStreamLine* streamer = vtkStreamLine::New();
333
#else
334
vtkStreamTracer* streamer = vtkStreamTracer::New();
335
#endif
336
337
338
vtkRungeKutta4* integrator = vtkRungeKutta4::New();
339
340
#if VTK_MAJOR_VERSION <= 5
341
streamer->SetInput(grid);
342
#else
343
streamer->SetInputData(grid);
344
#endif
345
346
if(lineSource) {
347
#if VTK_MAJOR_VERSION <= 5
348
streamer->SetSource(line->GetOutput());
349
#else
350
streamer->SetSourceConnection(line->GetOutputPort());
351
#endif
352
} else {
353
#if VTK_MAJOR_VERSION <= 5
354
streamer->SetSource(point->GetOutput());
355
#else
356
streamer->SetSourceConnection(point->GetOutputPort());
357
#endif
358
}
359
360
streamer->SetIntegrator(integrator);
361
362
#if VTK_MAJOR_VERSION <= 7
363
streamer->SetMaximumPropagationTime(propagationTime);
364
streamer->SetIntegrationStepLength(integStepLength);
365
366
if(forward && backward) {
367
streamer->SetIntegrationDirectionToIntegrateBothDirections();
368
} else if(forward) {
369
streamer->SetIntegrationDirectionToForward();
370
} else {
371
streamer->SetIntegrationDirectionToBackward();
372
}
373
streamer->SetStepLength(stepLength);
374
streamer->SetNumberOfThreads(threads);
375
376
#else
377
streamer->SetMaximumPropagation(propagationTime);
378
streamer->SetInitialIntegrationStep(stepLength);
379
380
if(forward && backward) {
381
streamer->SetIntegrationDirectionToBoth();
382
} else if(forward) {
383
streamer->SetIntegrationDirectionToForward();
384
} else {
385
streamer->SetIntegrationDirectionToBackward();
386
}
387
streamer->SetMaximumIntegrationStep(integStepLength);
388
streamer->SetMaximumNumberOfSteps(threads);
389
390
#endif
391
392
vtkRibbonFilter* ribbon = vtkRibbonFilter::New();
393
394
if(drawRibbon) {
395
double length = grid->GetLength();
396
ribbon->SetInputConnection(streamer->GetOutputPort());
397
ribbon->SetWidth(ribbonWidth * length / 1000.0);
398
ribbon->SetWidthFactor(5);
399
}
400
401
vtkPolyDataMapper* mapper = vtkPolyDataMapper::New();
402
403
mapper->ScalarVisibilityOn();
404
mapper->SetScalarRange(minVal, maxVal);
405
406
if(drawRibbon) {
407
mapper->SetInputConnection(ribbon->GetOutputPort());
408
} else {
409
mapper->SetInputConnection(streamer->GetOutputPort());
410
}
411
412
mapper->SetColorModeToMapScalars();
413
mapper->InterpolateScalarsBeforeMappingOn();
414
//mapper->SetLookupTable(vtkPost->GetCurrentLut());
415
mapper->SetLookupTable(vtkPost->GetLut("Streamline"));
416
if(ui.colorCombo->currentIndex() == 0 ){ // i.e. Null field
417
mapper->SetScalarRange(0, 1);
418
double h = nullColor.hueF();
419
double s = nullColor.saturationF();
420
double v = nullColor.valueF();
421
int nColor =128;
422
vtkLookupTable* nullLut = vtkLookupTable::New();
423
nullLut->SetHueRange(h, h);
424
nullLut->SetSaturationRange(s, s);
425
nullLut->SetValueRange(v, v);
426
nullLut->SetNumberOfColors(nColor);
427
nullLut->Build();
428
mapper->SetLookupTable(nullLut);
429
nullLut->Delete();
430
}
431
vtkPost->GetStreamLineActor()->SetMapper(mapper);
432
433
if(!drawRibbon)
434
vtkPost->GetStreamLineActor()->GetProperty()->SetLineWidth(lineWidth);
435
436
vtkPost->SetCurrentStreamLineName(colorName);
437
438
// Clean up:
439
//----------
440
line->Delete();
441
point->Delete();
442
vectorData->Delete();
443
vectorColor->Delete();
444
integrator->Delete();
445
streamer->Delete();
446
ribbon->Delete();
447
mapper->Delete();
448
}
449
450
// Public slots:
451
//---------------
452
QString StreamLine::GetFieldName()
453
{
454
return ui.vectorCombo->currentText();
455
}
456
457
QString StreamLine::GetColorName()
458
{
459
return ui.colorCombo->currentText();
460
}
461
462
bool StreamLine::SetFieldName(QString name)
463
{
464
for(int i = 0; i < ui.vectorCombo->count(); i++) {
465
if(ui.vectorCombo->itemText(i) == name) {
466
ui.vectorCombo->setCurrentIndex(i);
467
return true;
468
}
469
}
470
return false;
471
}
472
473
bool StreamLine::SetColorName(QString name)
474
{
475
for(int i = 0; i < ui.colorCombo->count(); i++) {
476
if(ui.colorCombo->itemText(i) == name) {
477
ui.colorCombo->setCurrentIndex(i);
478
return true;
479
}
480
}
481
return false;
482
}
483
484
void StreamLine::SetMaxTime(double f)
485
{
486
ui.propagationTime->setText(QString::number(f));
487
}
488
489
void StreamLine::SetStepLength(double f)
490
{
491
ui.stepLength->setText(QString::number(f));
492
}
493
494
void StreamLine::SetThreads(int n)
495
{
496
ui.threads->setValue(n);
497
}
498
499
void StreamLine::UseSurfaceMesh(bool b)
500
{
501
ui.useSurfaceGrid->setChecked(b);
502
}
503
504
void StreamLine::UseVolumeMesh(bool b)
505
{
506
ui.useVolumeGrid->setChecked(b);
507
}
508
509
void StreamLine::IntegrateForwards(bool b)
510
{
511
ui.forward->setChecked(b);
512
}
513
514
void StreamLine::IntegrateBackwards(bool b)
515
{
516
ui.backward->setChecked(b);
517
}
518
519
void StreamLine::SetIntegStepLength(double f)
520
{
521
ui.integStepLength->setText(QString::number(f));
522
}
523
524
void StreamLine::SetMinColorVal(double f)
525
{
526
ui.minVal->setText(QString::number(f));
527
}
528
529
void StreamLine::SetMaxColorVal(double f)
530
{
531
ui.maxVal->setText(QString::number(f));
532
}
533
534
void StreamLine::KeepColorLimits(bool b)
535
{
536
ui.keepLimits->setChecked(b);
537
}
538
539
void StreamLine::DrawLines(bool b)
540
{
541
ui.drawLines->setChecked(b);
542
}
543
544
void StreamLine::DrawRibbons(bool b)
545
{
546
ui.drawRibbon->setChecked(b);
547
}
548
549
void StreamLine::SetLineWidth(int n)
550
{
551
ui.lineWidth->setValue(n);
552
}
553
554
void StreamLine::SetRibbonWidth(int n)
555
{
556
ui.ribbonWidth->setValue(n);
557
}
558
559
void StreamLine::UseSphereSource(bool b)
560
{
561
ui.sphereSource->setChecked(b);
562
}
563
564
void StreamLine::UseLineSource(bool b)
565
{
566
ui.lineSource->setChecked(b);
567
}
568
569
void StreamLine::UsePointSource(bool b)
570
{
571
ui.pickSource->setChecked(b);
572
}
573
574
void StreamLine::SetSphereSourceX(double f)
575
{
576
ui.centerX->setText(QString::number(f));
577
}
578
579
void StreamLine::SetSphereSourceY(double f)
580
{
581
ui.centerY->setText(QString::number(f));
582
}
583
584
void StreamLine::SetSphereSourceZ(double f)
585
{
586
ui.centerZ->setText(QString::number(f));
587
}
588
589
void StreamLine::SetSphereSourceRadius(double f)
590
{
591
ui.radius->setText(QString::number(f));
592
}
593
594
void StreamLine::SetSphereSourcePoints(int n)
595
{
596
ui.points->setValue(n);
597
}
598
599
void StreamLine::SetLineSourceStartX(double f)
600
{
601
ui.startX->setText(QString::number(f));
602
}
603
604
void StreamLine::SetLineSourceStartY(double f)
605
{
606
ui.startY->setText(QString::number(f));
607
}
608
609
void StreamLine::SetLineSourceStartZ(double f)
610
{
611
ui.startZ->setText(QString::number(f));
612
}
613
614
void StreamLine::SetLineSourceEndX(double f)
615
{
616
ui.endX->setText(QString::number(f));
617
}
618
619
void StreamLine::SetLineSourceEndY(double f)
620
{
621
ui.endY->setText(QString::number(f));
622
}
623
624
void StreamLine::SetLineSourceEndZ(double f)
625
{
626
ui.endZ->setText(QString::number(f));
627
}
628
629
void StreamLine::SetLineSourcePoints(int n)
630
{
631
ui.lines->setValue(n);
632
}
633
634
void StreamLine::nullColorButtonClicked()
635
{
636
setNullColor(QColorDialog::getColor(nullColor));
637
}
638
639
void StreamLine::setNullColor(QColor color){
640
if(!color.isValid()) return;
641
642
nullColor = color;
643
644
QPalette plt(ui.nullColorLabel->palette());
645
plt.setColor(QPalette::WindowText, nullColor);
646
ui.nullColorLabel->setPalette(plt);
647
}
648
649