OGS
VtkVisPipeline.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4// ** INCLUDES **
5#include "VtkVisPipeline.h"
6
7#include <vtkAlgorithm.h>
8#include <vtkCamera.h>
9#include <vtkCellData.h>
10#include <vtkFieldData.h>
11#include <vtkGenericDataObjectReader.h>
12#include <vtkImageActor.h>
13#include <vtkImageReader2.h>
14#include <vtkLight.h>
15#include <vtkPointData.h>
16#include <vtkPointSet.h>
17#include <vtkProp3D.h>
18#include <vtkRenderer.h>
19#include <vtkSmartPointer.h>
20#include <vtkTransformFilter.h>
21#include <vtkXMLImageDataReader.h>
22#include <vtkXMLPolyDataReader.h>
23#include <vtkXMLRectilinearGridReader.h>
24#include <vtkXMLStructuredGridReader.h>
25#include <vtkXMLUnstructuredGridReader.h>
26
27#include <QColor>
28#ifndef NDEBUG
29#include <QElapsedTimer>
30#endif // NDEBUG
31#include <QFileInfo>
32#include <QSettings>
33#include <QString>
34
35#include "Base/TreeModel.h"
36#include "BaseLib/Logging.h"
38#include "DataView/MeshItem.h"
39#include "DataView/MeshModel.h"
42#include "MeshLib/Mesh.h"
48#include "VtkFilterFactory.h"
49#include "VtkVisImageItem.h"
50#include "VtkVisPipelineItem.h"
51#include "VtkVisPointSetItem.h"
52
53VtkVisPipeline::VtkVisPipeline(vtkRenderer* renderer, QObject* parent /*= 0*/)
54 : TreeModel(parent), _renderer(renderer)
55{
56 QList<QVariant> rootData;
57 rootData << "Object name"
58 << "Visible";
59 delete _rootItem;
60 _rootItem = new TreeItem(rootData, nullptr);
61
62 QSettings settings;
63 QVariant backgroundColorVariant = settings.value("VtkBackgroundColor");
64 if (backgroundColorVariant != QVariant())
65 {
66 this->setBGColor(backgroundColorVariant.value<QColor>());
67 }
68
70 settings.value("resetViewOnLoad", true).toBool();
71}
72
73bool VtkVisPipeline::setData(const QModelIndex& index, const QVariant& value,
74 int role /* = Qt::EditRole */)
75{
77
78 return TreeModel::setData(index, value, role);
79}
80
82{
83 double lightPos[3];
84 for (auto& light : _lights)
85 {
86 light->GetPosition(lightPos);
87 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
88 pos[2] == lightPos[2])
89 {
90 return;
91 }
92 }
93 vtkLight* l = vtkLight::New();
94 l->SetPosition(pos[0], pos[1], pos[2]);
95 _renderer->AddLight(l);
96 _lights.push_back(l);
97}
98
99vtkLight* VtkVisPipeline::getLight(const GeoLib::Point& pos) const
100{
101 double lightPos[3];
102 for (auto light : _lights)
103 {
104 light->GetPosition(lightPos);
105 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
106 pos[2] == lightPos[2])
107 {
108 return light;
109 }
110 }
111 return nullptr;
112}
113
115{
116 double lightPos[3];
117 for (auto it = _lights.begin(); it != _lights.end(); ++it)
118 {
119 (*it)->GetPosition(lightPos);
120 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
121 pos[2] == lightPos[2])
122 {
123 _renderer->RemoveLight(*it);
124 (*it)->Delete();
125 _lights.erase(it);
126 return;
127 }
128 }
129}
130
132{
133 double* color = _renderer->GetBackground();
134 QColor c(static_cast<int>(color[0] * 255),
135 static_cast<int>(color[1] * 255),
136 static_cast<int>(color[2] * 255));
137 return c;
138}
139
140void VtkVisPipeline::setBGColor(const QColor& color)
141{
142 QSettings settings;
143 settings.setValue("VtkBackgroundColor", color);
144 _renderer->SetBackground(color.redF(), color.greenF(), color.blueF());
145}
146
147QModelIndex VtkVisPipeline::getIndex(vtkProp3D* actor)
148{
149 return _actorMap.value(actor, QModelIndex());
150}
151
152Qt::ItemFlags VtkVisPipeline::flags(const QModelIndex& index) const
153{
154 Qt::ItemFlags defaultFlags = Qt::ItemIsEnabled | Qt::ItemIsSelectable;
155
156 if (!index.isValid())
157 {
158 return Qt::ItemIsEnabled;
159 }
160
161 // if (index.column() == 1)
162 // defaultFlags |= Qt::ItemIsEditable;
163
164 return defaultFlags;
165}
166
167void VtkVisPipeline::loadFromFile(QString filename)
168{
169#ifndef NDEBUG
170 QElapsedTimer myTimer;
171 myTimer.start();
172 INFO("VTK Read: {:s}.", filename.toStdString());
173#endif
174
175 if (filename.size() > 0)
176 {
177 vtkSmartPointer<vtkXMLDataReader> reader;
178 if (filename.endsWith("vti"))
179 {
180 reader = vtkSmartPointer<vtkXMLImageDataReader>::New();
181 }
182 else if (filename.endsWith("vtr"))
183 {
184 reader = vtkSmartPointer<vtkXMLRectilinearGridReader>::New();
185 }
186 else if (filename.endsWith("vts"))
187 {
188 reader = vtkSmartPointer<vtkXMLStructuredGridReader>::New();
189 }
190 else if (filename.endsWith("vtp"))
191 {
192 reader = vtkSmartPointer<vtkXMLPolyDataReader>::New();
193 }
194 else if (filename.endsWith("vtu"))
195 {
196 reader = vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
197 }
198 else if (filename.endsWith("vtk"))
199 {
200 vtkGenericDataObjectReader* oldStyleReader =
201 vtkGenericDataObjectReader::New();
202 oldStyleReader->SetFileName(filename.toStdString().c_str());
203 oldStyleReader->ReadAllFieldsOn();
204 oldStyleReader->ReadAllScalarsOn();
205 oldStyleReader->Update();
206 vtkDataSet* dataSet =
207 vtkDataSet::SafeDownCast(oldStyleReader->GetOutput());
208 if (dataSet)
209 {
210 this->listArrays(dataSet);
211 addPipelineItem(oldStyleReader);
212 }
213 else
214 ERR("VtkVisPipeline::loadFromFile(): not a valid vtkDataSet.");
215 return;
216 }
217 else
218 {
219 return;
220 }
221
222 reader->SetFileName(filename.toStdString().c_str());
223 // TODO: insert ReadAllScalarsOn()-equivalent for xml-file-reader here,
224 // otherwise arrays are not available in GUI!
225 reader->Update();
226 // std::cout << "#cell scalars: " << reader->GetNumberOfCellArrays() <<
227 // std::endl; std::cout << "#point scalars: " <<
228 // reader->GetNumberOfPointArrays() << std::endl;
229
230 vtkSmartPointer<vtkDataSet> dataSet = reader->GetOutputAsDataSet();
231 if (dataSet)
232 {
233 this->listArrays(dataSet);
234 addPipelineItem(reader);
235 }
236 else
237 ERR("VtkVisPipeline::loadFromFile(): not a valid vtkDataSet.");
238 }
239
240#ifndef NDEBUG
241 INFO("{:d} ms", myTimer.elapsed());
242#endif
243}
244
246{
247 // iterate over all source items
248 for (int i = 0; i < _rootItem->childCount(); ++i)
249 {
250 auto* item = static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
251 item->setScale(1.0, 1.0, factor);
252
253 // recursively set on all child items
254 item->setScaleOnChildren(1.0, 1.0, 1.0);
255 }
256
258}
259
261{
262 // iterate over all source items
263 for (int i = 0; i < _rootItem->childCount(); ++i)
264 {
265 auto* item = static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
266 item->setBackfaceCulling(enable);
267
268 // recursively set on all child items
269 item->setBackfaceCullingOnChildren(enable);
270 }
271
273}
274
276 const std::string& name,
277 GeoLib::GEOTYPE type)
278{
279 addPipelineItem(model->vtkSource(name, type));
280}
281
283 const std::string& name)
284{
285 addPipelineItem(model->vtkSource(name));
286}
287
288void VtkVisPipeline::addPipelineItem(MeshModel* model, const QModelIndex& idx)
289{
290 addPipelineItem(static_cast<MeshItem*>(model->getItem(idx))->vtkSource());
291}
292
294 const QModelIndex& parent)
295{
296 beginResetModel();
297
298 item->Initialize(_renderer);
299 TreeItem* parentItem = item->parentItem();
300 parentItem->appendChild(item);
301
302 if (!parent.isValid()) // Set global superelevation on source objects
303 {
304 QSettings settings;
305 if (dynamic_cast<vtkImageAlgorithm*>(item->algorithm()) == nullptr)
306 { // if not an image
307 item->setScale(
308 1.0, 1.0,
309 settings.value("globalSuperelevation", 1.0).toDouble());
310 }
311 }
312
313 int parentChildCount = parentItem->childCount();
314 QModelIndex newIndex = index(parentChildCount - 1, 0, parent);
315
316 if (_resetCameraOnAddOrRemove || _rootItem->childCount() == 1)
317 {
318 _renderer->ResetCamera(_renderer->ComputeVisiblePropBounds());
319 }
320 _actorMap.insert(item->actor(), newIndex);
321
322 // Do not interpolate images
323 if (dynamic_cast<vtkImageAlgorithm*>(item->algorithm()))
324 {
325 static_cast<vtkImageActor*>(item->actor())->InterpolateOff();
326 }
327
328 endResetModel();
330 emit itemSelected(newIndex);
331
332 return newIndex;
333}
334
336 vtkAlgorithm* source, QModelIndex parent /* = QModelindex() */)
337{
338 std::string itemName;
339
340 if (!parent.isValid()) // if source object
341 {
342 auto* old_reader = dynamic_cast<vtkGenericDataObjectReader*>(source);
343 auto* new_reader = dynamic_cast<vtkXMLReader*>(source);
344 auto* image_reader = dynamic_cast<vtkImageReader2*>(source);
345 auto* props = dynamic_cast<VtkAlgorithmProperties*>(source);
346 auto* meshSource = dynamic_cast<MeshLib::VtkMappedMeshSource*>(source);
347 if (old_reader)
348 {
349 itemName = old_reader->GetFileName();
350 }
351 else if (new_reader)
352 {
353 itemName = new_reader->GetFileName();
354 }
355 else if (image_reader)
356 {
357 itemName = image_reader->GetFileName();
358 }
359 else if (props)
360 {
361 itemName = props->GetName().toStdString();
362 }
363 else if (meshSource)
364 {
365 itemName = meshSource->GetMesh()->getName();
366 }
367 }
368
369 if (itemName.length() == 0)
370 {
371 itemName = source->GetClassName();
372 }
373
374 QList<QVariant> itemData;
375 itemData << QString::fromStdString(itemName) << true;
376
377 VtkVisPipelineItem* item(nullptr);
378 if (dynamic_cast<vtkImageAlgorithm*>(source))
379 {
380 item = new VtkVisImageItem(source, getItem(parent), itemData);
381 }
382 else
383 {
384 item = new VtkVisPointSetItem(source, getItem(parent), itemData);
385 }
386 return this->addPipelineItem(item, parent);
387}
388
390 const std::string& name,
391 GeoLib::GEOTYPE type)
392{
393 for (int i = 0; i < _rootItem->childCount(); i++)
394 {
395 VtkVisPipelineItem* item =
396 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
397 if (item->algorithm() == model->vtkSource(name, type))
398 {
400 return;
401 }
402 }
403}
404
406 const std::string& name)
407{
408 for (int i = 0; i < _rootItem->childCount(); i++)
409 {
410 VtkVisPipelineItem* item =
411 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
412 if (item->algorithm() == model->vtkSource(name))
413 {
415 return;
416 }
417 }
418}
419
420void VtkVisPipeline::removeSourceItem(MeshModel* model, const QModelIndex& idx)
421{
422 auto* sItem = static_cast<MeshItem*>(model->getItem(idx));
423
424 for (int i = 0; i < _rootItem->childCount(); i++)
425 {
426 VtkVisPipelineItem* item =
427 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
428 if (item->algorithm() == sItem->vtkSource())
429 {
431 return;
432 }
433 }
434}
435
437{
438 if (!index.isValid())
439 {
440 return;
441 }
442
443 QMap<vtkProp3D*, QModelIndex>::iterator it = _actorMap.begin();
444 while (it != _actorMap.end())
445 {
446 QModelIndex itIndex = it.value();
447 if (itIndex == index)
448 {
449 _actorMap.erase(it);
450 break;
451 }
452 ++it;
453 }
454
455 // TreeItem* item = getItem(index);
456 removeRows(index.row(), 1, index.parent());
457
459 {
460 _renderer->ResetCamera(_renderer->ComputeVisiblePropBounds());
461 }
463}
464
465void VtkVisPipeline::listArrays(vtkDataSet* dataSet)
466{
467 if (dataSet)
468 {
469 vtkPointData* pointData = dataSet->GetPointData();
470 INFO(" #point data arrays: {:d}", pointData->GetNumberOfArrays());
471 for (int i = 0; i < pointData->GetNumberOfArrays(); i++)
472 INFO(" Name: {:s}", pointData->GetArrayName(i));
473
474 vtkCellData* cellData = dataSet->GetCellData();
475 INFO(" #cell data arrays: {:d}", cellData->GetNumberOfArrays());
476 for (int i = 0; i < cellData->GetNumberOfArrays(); i++)
477 INFO(" Name: {:s}", cellData->GetArrayName(i));
478 }
479 else
480 ERR("VtkVisPipeline::listArrays(): not a valid vtkDataSet.");
481}
482
485 std::vector<double> const& quality)
486{
487 if (!source || quality.empty())
488 {
489 return;
490 }
491
492 int const nSources = this->_rootItem->childCount();
493 for (int i = 0; i < nSources; i++)
494 {
495 auto* parentItem =
496 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
497 if (parentItem->algorithm() != source)
498 {
499 continue;
500 }
501
502 QList<QVariant> itemData;
503 itemData << "MeshQuality: " +
504 QString::fromStdString(MeshQualityType2String(t))
505 << true;
506
508 "VtkCompositeElementSelectionFilter",
509 parentItem->transformFilter());
511 {
512 auto const range(
513 std::minmax_element(quality.cbegin(), quality.cend()));
514 static_cast<VtkCompositeElementSelectionFilter*>(filter)->setRange(
515 *range.first, *range.second);
516 }
517 static_cast<VtkCompositeElementSelectionFilter*>(filter)
518 ->setSelectionArray("Selection", quality);
519 VtkVisPointSetItem* item =
520 new VtkVisPointSetItem(filter, parentItem, itemData);
521 this->addPipelineItem(item, this->createIndex(i, 0, item));
522 break;
523 }
524}
525
526void VtkVisPipeline::highlightGeoObject(const vtkPolyDataAlgorithm* source,
527 int index)
528{
530 int nSources = this->_rootItem->childCount();
531 for (int i = 0; i < nSources; i++)
532 {
533 auto* parentItem =
534 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
535 if (parentItem->algorithm() == source)
536 {
537 QList<QVariant> itemData;
538 itemData << "Selected GeoObject" << true;
539
540 VtkCompositeFilter* filter =
542 "VtkCompositeGeoObjectFilter",
543 parentItem->transformFilter());
544 static_cast<VtkCompositeGeoObjectFilter*>(filter)->SetIndex(index);
545 VtkVisPointSetItem* item =
546 new VtkVisPointSetItem(filter, parentItem, itemData);
547 QModelIndex parent_index =
548 static_cast<TreeModel*>(this)->index(i, 0, QModelIndex());
549 _highlighted_geo_index = this->addPipelineItem(item, parent_index);
550 break;
551 }
552 }
553}
554
556{
557 if (_highlighted_geo_index != QModelIndex())
558 {
560 _highlighted_geo_index = QModelIndex();
561 }
562}
563
565 vtkUnstructuredGridAlgorithm const* const source, unsigned index,
566 bool is_element)
567{
568 int nSources = this->_rootItem->childCount();
569 for (int i = 0; i < nSources; i++)
570 {
571 auto* parentItem =
572 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
573 if (parentItem->algorithm() == source)
574 {
575 QList<QVariant> itemData;
576 itemData << "Selected Mesh Component" << true;
577 QList<QVariant> selected_index;
578 selected_index << index << index;
579
580 VtkCompositeFilter* filter(nullptr);
581 if (is_element)
582 {
584 "VtkCompositeElementSelectionFilter",
585 parentItem->transformFilter());
586 static_cast<VtkCompositeElementSelectionFilter*>(filter)
587 ->setSelectionArray("vtkIdFilter_Ids");
588 static_cast<VtkCompositeElementSelectionFilter*>(filter)
589 ->SetUserVectorProperty("Threshold Between",
590 selected_index);
591 }
592 else
593 {
595 "VtkCompositeNodeSelectionFilter",
596 parentItem->transformFilter());
597 std::vector<unsigned> indices(1);
598 indices[0] = index;
599 static_cast<VtkCompositeNodeSelectionFilter*>(filter)
600 ->setSelectionArray(indices);
601 }
602 VtkVisPointSetItem* item =
603 new VtkVisPointSetItem(filter, parentItem, itemData);
604 QModelIndex parent_index =
605 static_cast<TreeModel*>(this)->index(i, 0, QModelIndex());
607 this->addPipelineItem(item, parent_index);
608 break;
609 }
610 }
611}
612
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
A model for the GeoTreeView implementing a tree as a double-linked list.
vtkPolyDataAlgorithm * vtkSource(const std::string &name, GeoLib::GEOTYPE type) const
Returns the vtk-object indicated by type of the geometry indicated by name.
A TreeItem containing a mesh and the associated vtk object used in the Mesh Model.
Definition MeshItem.h:19
MeshLib::VtkMappedMeshSource * vtkSource() const
Returns the VTK object.
Definition MeshItem.h:33
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
A model for the StationTreeView implementing a tree as a double-linked list.
vtkPolyDataAlgorithm * vtkSource(const std::string &name) const
Objects nodes for the TreeModel.
Definition TreeItem.h:17
virtual int childCount() const
Definition TreeItem.cpp:54
void appendChild(TreeItem *item)
Definition TreeItem.cpp:31
TreeItem * parentItem() const
Definition TreeItem.cpp:104
bool removeRows(int position, int count, const QModelIndex &parent) override
TreeItem * getItem(const QModelIndex &index) const
QModelIndex index(int row, int column, const QModelIndex &parent=QModelIndex()) const override
Definition TreeModel.cpp:39
QModelIndex parent(const QModelIndex &index) const override
Definition TreeModel.cpp:70
TreeModel(QObject *parent=nullptr)
Definition TreeModel.cpp:15
TreeItem * _rootItem
Definition TreeModel.h:47
bool setData(const QModelIndex &index, const QVariant &value, int role) override
Contains properties for the visualization of objects as VtkVisPipelineItems.
This filter selects/thresholds elements based on the selected array.
Is used to combine several filter in one VtkVisPipelineItem. You can use vtk filter and custom filter...
Highlights a single GeoObject.
This filter displays the points/nodes given in the index field as spheres.
static VtkCompositeFilter * CreateCompositeFilter(QString type, vtkAlgorithm *inputAlgorithm)
Creates a composite filter by name.
An item in the VtkVisPipeline containing an image to be visualized.
An item in the VtkVisPipeline containing a graphic object to be visualized.
virtual void setScale(double x, double y, double z) const
Scales the data in visualisation-space. This function is empty and needs to be implemented by derived...
vtkProp3D * actor() const
Returns the actor as vtkProp3D.
vtkAlgorithm * algorithm() const
Returns the algorithm object.
virtual void setBackfaceCulling(bool enable) const
Enables / disables backface culling.
virtual void Initialize(vtkRenderer *renderer)=0
Initializes vtkMapper and vtkActor necessary for visualization of the item and sets the item's proper...
void setGlobalBackfaceCulling(bool enable) const
Enables / disables backface culling on all actors.
QModelIndex _highlighted_mesh_component
void itemSelected(const QModelIndex &) const
bool setData(const QModelIndex &index, const QVariant &value, int role=Qt::EditRole) override
Emits vtkVisPipelineChanged() and calls base class method.
std::list< vtkLight * > _lights
void addLight(const GeoLib::Point &pos)
Adds a light to the scene at the given coordinates.
void removeLight(const GeoLib::Point &pos)
Removes a light at the given coordinates (if possible).
void vtkVisPipelineChanged() const
Is emitted when a pipeline item was added or removed.
void addPipelineItem(MeshModel *model, const QModelIndex &idx)
Adds the given Model to the pipeline.
void removePipelineItem(QModelIndex index)
Removes the vtkAlgorithm at the given QModelIndex (and all attached vtkAlgorithms) from the pipeline.
VtkVisPipeline(vtkRenderer *renderer, QObject *parent=nullptr)
void removeHighlightedGeoObject()
Removes the currently highlighted geometry-object.
void setBGColor(const QColor &color)
Sets the background-colour of the scene.
void setGlobalSuperelevation(double factor) const
Sets a global superelevation factor on all source items and resets the factor on other items to 1.
void removeHighlightedMeshComponent()
Removes the currently highlighted mesh component.
Qt::ItemFlags flags(const QModelIndex &index) const override
vtkLight * getLight(const GeoLib::Point &pos) const
Returns a light (or nullptr) for the given coordinates.
void loadFromFile(QString filename)
Loads a vtk object from the given file and adds it to the pipeline.
void highlightGeoObject(const vtkPolyDataAlgorithm *source, int index)
Applies a VtkCompositeGeoObjectFilter to add a specific index of the given geometry-source to the pip...
void highlightMeshComponent(vtkUnstructuredGridAlgorithm const *const source, unsigned index, bool is_element)
Applies a VtkCompositeSelectionFilter to add a specific component of the given mesh-source to the pip...
void showMeshElementQuality(MeshLib::VtkMappedMeshSource *source, MeshLib::MeshQualityType t, std::vector< double > const &quality)
Checks the quality of mesh elements and adds a filter to highlight deformed elements.
void removeSourceItem(MeshModel *model, const QModelIndex &idx)
Removes the given Model (and all attached vtkAlgorithms) from the pipeline.
QModelIndex getIndex(vtkProp3D *actor)
Returns the QModelIndex of VtkVisPipelineItem which actor is the given one.
bool _resetCameraOnAddOrRemove
QMap< vtkProp3D *, QModelIndex > _actorMap
QModelIndex _highlighted_geo_index
QColor getBGColor() const
Returns the background-colour of the scene.
void listArrays(vtkDataSet *dataSet)
vtkRenderer * _renderer
An item in the VtkVisPipeline containing a point set object to be visualized.
GEOTYPE
Definition GeoType.h:12
MeshQualityType
Describes a mesh quality metric.
Definition MeshEnums.h:80