OGS
VtkVisPipeline.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
16#include "VtkVisPipeline.h"
17
18#include <vtkAlgorithm.h>
19#include <vtkCamera.h>
20#include <vtkCellData.h>
21#include <vtkFieldData.h>
22#include <vtkGenericDataObjectReader.h>
23#include <vtkImageActor.h>
24#include <vtkImageReader2.h>
25#include <vtkLight.h>
26#include <vtkPointData.h>
27#include <vtkPointSet.h>
28#include <vtkProp3D.h>
29#include <vtkRenderer.h>
30#include <vtkSmartPointer.h>
31#include <vtkTransformFilter.h>
32#include <vtkXMLImageDataReader.h>
33#include <vtkXMLPolyDataReader.h>
34#include <vtkXMLRectilinearGridReader.h>
35#include <vtkXMLStructuredGridReader.h>
36#include <vtkXMLUnstructuredGridReader.h>
37
38#include <QColor>
39#ifndef NDEBUG
40#include <QElapsedTimer>
41#endif // NDEBUG
42#include <QFileInfo>
43#include <QSettings>
44#include <QString>
45
46#include "Base/TreeModel.h"
47#include "BaseLib/Logging.h"
49#include "DataView/MeshItem.h"
50#include "DataView/MeshModel.h"
53#include "MeshLib/Mesh.h"
59#include "VtkFilterFactory.h"
60#include "VtkVisImageItem.h"
61#include "VtkVisPipelineItem.h"
62#include "VtkVisPointSetItem.h"
63
64VtkVisPipeline::VtkVisPipeline(vtkRenderer* renderer, QObject* parent /*= 0*/)
65 : TreeModel(parent), _renderer(renderer)
66{
67 QList<QVariant> rootData;
68 rootData << "Object name"
69 << "Visible";
70 delete _rootItem;
71 _rootItem = new TreeItem(rootData, nullptr);
72
73 QSettings settings;
74 QVariant backgroundColorVariant = settings.value("VtkBackgroundColor");
75 if (backgroundColorVariant != QVariant())
76 {
77 this->setBGColor(backgroundColorVariant.value<QColor>());
78 }
79
81 settings.value("resetViewOnLoad", true).toBool();
82}
83
84bool VtkVisPipeline::setData(const QModelIndex& index, const QVariant& value,
85 int role /* = Qt::EditRole */)
86{
88
89 return TreeModel::setData(index, value, role);
90}
91
93{
94 double lightPos[3];
95 for (auto& light : _lights)
96 {
97 light->GetPosition(lightPos);
98 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
99 pos[2] == lightPos[2])
100 {
101 return;
102 }
103 }
104 vtkLight* l = vtkLight::New();
105 l->SetPosition(pos[0], pos[1], pos[2]);
106 _renderer->AddLight(l);
107 _lights.push_back(l);
108}
109
110vtkLight* VtkVisPipeline::getLight(const GeoLib::Point& pos) const
111{
112 double lightPos[3];
113 for (auto light : _lights)
114 {
115 light->GetPosition(lightPos);
116 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
117 pos[2] == lightPos[2])
118 {
119 return light;
120 }
121 }
122 return nullptr;
123}
124
126{
127 double lightPos[3];
128 for (auto it = _lights.begin(); it != _lights.end(); ++it)
129 {
130 (*it)->GetPosition(lightPos);
131 if (pos[0] == lightPos[0] && pos[1] == lightPos[1] &&
132 pos[2] == lightPos[2])
133 {
134 _renderer->RemoveLight(*it);
135 (*it)->Delete();
136 _lights.erase(it);
137 return;
138 }
139 }
140}
141
143{
144 double* color = _renderer->GetBackground();
145 QColor c(static_cast<int>(color[0] * 255),
146 static_cast<int>(color[1] * 255),
147 static_cast<int>(color[2] * 255));
148 return c;
149}
150
151void VtkVisPipeline::setBGColor(const QColor& color)
152{
153 QSettings settings;
154 settings.setValue("VtkBackgroundColor", color);
155 _renderer->SetBackground(color.redF(), color.greenF(), color.blueF());
156}
157
158QModelIndex VtkVisPipeline::getIndex(vtkProp3D* actor)
159{
160 return _actorMap.value(actor, QModelIndex());
161}
162
163Qt::ItemFlags VtkVisPipeline::flags(const QModelIndex& index) const
164{
165 Qt::ItemFlags defaultFlags = Qt::ItemIsEnabled | Qt::ItemIsSelectable;
166
167 if (!index.isValid())
168 {
169 return Qt::ItemIsEnabled;
170 }
171
172 // if (index.column() == 1)
173 // defaultFlags |= Qt::ItemIsEditable;
174
175 return defaultFlags;
176}
177
178void VtkVisPipeline::loadFromFile(QString filename)
179{
180#ifndef NDEBUG
181 QElapsedTimer myTimer;
182 myTimer.start();
183 INFO("VTK Read: {:s}.", filename.toStdString());
184#endif
185
186 if (filename.size() > 0)
187 {
188 vtkSmartPointer<vtkXMLDataReader> reader;
189 if (filename.endsWith("vti"))
190 {
191 reader = vtkSmartPointer<vtkXMLImageDataReader>::New();
192 }
193 else if (filename.endsWith("vtr"))
194 {
195 reader = vtkSmartPointer<vtkXMLRectilinearGridReader>::New();
196 }
197 else if (filename.endsWith("vts"))
198 {
199 reader = vtkSmartPointer<vtkXMLStructuredGridReader>::New();
200 }
201 else if (filename.endsWith("vtp"))
202 {
203 reader = vtkSmartPointer<vtkXMLPolyDataReader>::New();
204 }
205 else if (filename.endsWith("vtu"))
206 {
207 reader = vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
208 }
209 else if (filename.endsWith("vtk"))
210 {
211 vtkGenericDataObjectReader* oldStyleReader =
212 vtkGenericDataObjectReader::New();
213 oldStyleReader->SetFileName(filename.toStdString().c_str());
214 oldStyleReader->ReadAllFieldsOn();
215 oldStyleReader->ReadAllScalarsOn();
216 oldStyleReader->Update();
217 vtkDataSet* dataSet =
218 vtkDataSet::SafeDownCast(oldStyleReader->GetOutput());
219 if (dataSet)
220 {
221 this->listArrays(dataSet);
222 addPipelineItem(oldStyleReader);
223 }
224 else
225 ERR("VtkVisPipeline::loadFromFile(): not a valid vtkDataSet.");
226 return;
227 }
228 else
229 {
230 return;
231 }
232
233 reader->SetFileName(filename.toStdString().c_str());
234 // TODO: insert ReadAllScalarsOn()-equivalent for xml-file-reader here,
235 // otherwise arrays are not available in GUI!
236 reader->Update();
237 // std::cout << "#cell scalars: " << reader->GetNumberOfCellArrays() <<
238 // std::endl; std::cout << "#point scalars: " <<
239 // reader->GetNumberOfPointArrays() << std::endl;
240
241 vtkSmartPointer<vtkDataSet> dataSet = reader->GetOutputAsDataSet();
242 if (dataSet)
243 {
244 this->listArrays(dataSet);
245 addPipelineItem(reader);
246 }
247 else
248 ERR("VtkVisPipeline::loadFromFile(): not a valid vtkDataSet.");
249 }
250
251#ifndef NDEBUG
252 INFO("{:d} ms", myTimer.elapsed());
253#endif
254}
255
257{
258 // iterate over all source items
259 for (int i = 0; i < _rootItem->childCount(); ++i)
260 {
261 auto* item = static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
262 item->setScale(1.0, 1.0, factor);
263
264 // recursively set on all child items
265 item->setScaleOnChildren(1.0, 1.0, 1.0);
266 }
267
269}
270
272{
273 // iterate over all source items
274 for (int i = 0; i < _rootItem->childCount(); ++i)
275 {
276 auto* item = static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
277 item->setBackfaceCulling(enable);
278
279 // recursively set on all child items
280 item->setBackfaceCullingOnChildren(enable);
281 }
282
284}
285
287 const std::string& name,
288 GeoLib::GEOTYPE type)
289{
290 addPipelineItem(model->vtkSource(name, type));
291}
292
294 const std::string& name)
295{
296 addPipelineItem(model->vtkSource(name));
297}
298
299void VtkVisPipeline::addPipelineItem(MeshModel* model, const QModelIndex& idx)
300{
301 addPipelineItem(static_cast<MeshItem*>(model->getItem(idx))->vtkSource());
302}
303
305 const QModelIndex& parent)
306{
307 beginResetModel();
308
309 item->Initialize(_renderer);
310 TreeItem* parentItem = item->parentItem();
311 parentItem->appendChild(item);
312
313 if (!parent.isValid()) // Set global superelevation on source objects
314 {
315 QSettings settings;
316 if (dynamic_cast<vtkImageAlgorithm*>(item->algorithm()) == nullptr)
317 { // if not an image
318 item->setScale(
319 1.0, 1.0,
320 settings.value("globalSuperelevation", 1.0).toDouble());
321 }
322 }
323
324 int parentChildCount = parentItem->childCount();
325 QModelIndex newIndex = index(parentChildCount - 1, 0, parent);
326
328 {
329 _renderer->ResetCamera(_renderer->ComputeVisiblePropBounds());
330 }
331 _actorMap.insert(item->actor(), newIndex);
332
333 // Do not interpolate images
334 if (dynamic_cast<vtkImageAlgorithm*>(item->algorithm()))
335 {
336 static_cast<vtkImageActor*>(item->actor())->InterpolateOff();
337 }
338
339 endResetModel();
341 emit itemSelected(newIndex);
342
343 return newIndex;
344}
345
347 vtkAlgorithm* source, QModelIndex parent /* = QModelindex() */)
348{
349 std::string itemName;
350
351 if (!parent.isValid()) // if source object
352 {
353 auto* old_reader = dynamic_cast<vtkGenericDataObjectReader*>(source);
354 auto* new_reader = dynamic_cast<vtkXMLReader*>(source);
355 auto* image_reader = dynamic_cast<vtkImageReader2*>(source);
356 auto* props = dynamic_cast<VtkAlgorithmProperties*>(source);
357 auto* meshSource = dynamic_cast<MeshLib::VtkMappedMeshSource*>(source);
358 if (old_reader)
359 {
360 itemName = old_reader->GetFileName();
361 }
362 else if (new_reader)
363 {
364 itemName = new_reader->GetFileName();
365 }
366 else if (image_reader)
367 {
368 itemName = image_reader->GetFileName();
369 }
370 else if (props)
371 {
372 itemName = props->GetName().toStdString();
373 }
374 else if (meshSource)
375 {
376 itemName = meshSource->GetMesh()->getName();
377 }
378 }
379
380 if (itemName.length() == 0)
381 {
382 itemName = source->GetClassName();
383 }
384
385 QList<QVariant> itemData;
386 itemData << QString::fromStdString(itemName) << true;
387
388 VtkVisPipelineItem* item(nullptr);
389 if (dynamic_cast<vtkImageAlgorithm*>(source))
390 {
391 item = new VtkVisImageItem(source, getItem(parent), itemData);
392 }
393 else
394 {
395 item = new VtkVisPointSetItem(source, getItem(parent), itemData);
396 }
397 return this->addPipelineItem(item, parent);
398}
399
401 const std::string& name,
402 GeoLib::GEOTYPE type)
403{
404 for (int i = 0; i < _rootItem->childCount(); i++)
405 {
406 VtkVisPipelineItem* item =
407 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
408 if (item->algorithm() == model->vtkSource(name, type))
409 {
411 return;
412 }
413 }
414}
415
417 const std::string& name)
418{
419 for (int i = 0; i < _rootItem->childCount(); i++)
420 {
421 VtkVisPipelineItem* item =
422 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
423 if (item->algorithm() == model->vtkSource(name))
424 {
426 return;
427 }
428 }
429}
430
431void VtkVisPipeline::removeSourceItem(MeshModel* model, const QModelIndex& idx)
432{
433 auto* sItem = static_cast<MeshItem*>(model->getItem(idx));
434
435 for (int i = 0; i < _rootItem->childCount(); i++)
436 {
437 VtkVisPipelineItem* item =
438 static_cast<VtkVisPipelineItem*>(getItem(index(i, 0)));
439 if (item->algorithm() == sItem->vtkSource())
440 {
442 return;
443 }
444 }
445}
446
448{
449 if (!index.isValid())
450 {
451 return;
452 }
453
454 QMap<vtkProp3D*, QModelIndex>::iterator it = _actorMap.begin();
455 while (it != _actorMap.end())
456 {
457 QModelIndex itIndex = it.value();
458 if (itIndex == index)
459 {
460 _actorMap.erase(it);
461 break;
462 }
463 ++it;
464 }
465
466 // TreeItem* item = getItem(index);
467 removeRows(index.row(), 1, index.parent());
468
470 {
471 _renderer->ResetCamera(_renderer->ComputeVisiblePropBounds());
472 }
474}
475
476void VtkVisPipeline::listArrays(vtkDataSet* dataSet)
477{
478 if (dataSet)
479 {
480 vtkPointData* pointData = dataSet->GetPointData();
481 INFO(" #point data arrays: {:d}", pointData->GetNumberOfArrays());
482 for (int i = 0; i < pointData->GetNumberOfArrays(); i++)
483 INFO(" Name: {:s}", pointData->GetArrayName(i));
484
485 vtkCellData* cellData = dataSet->GetCellData();
486 INFO(" #cell data arrays: {:d}", cellData->GetNumberOfArrays());
487 for (int i = 0; i < cellData->GetNumberOfArrays(); i++)
488 INFO(" Name: {:s}", cellData->GetArrayName(i));
489 }
490 else
491 ERR("VtkVisPipeline::listArrays(): not a valid vtkDataSet.");
492}
493
496 std::vector<double> const& quality)
497{
498 if (!source || quality.empty())
499 {
500 return;
501 }
502
503 int const nSources = this->_rootItem->childCount();
504 for (int i = 0; i < nSources; i++)
505 {
506 auto* parentItem =
507 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
508 if (parentItem->algorithm() != source)
509 {
510 continue;
511 }
512
513 QList<QVariant> itemData;
514 itemData << "MeshQuality: " +
515 QString::fromStdString(MeshQualityType2String(t))
516 << true;
517
519 "VtkCompositeElementSelectionFilter",
520 parentItem->transformFilter());
522 {
523 auto const range(
524 std::minmax_element(quality.cbegin(), quality.cend()));
525 static_cast<VtkCompositeElementSelectionFilter*>(filter)->setRange(
526 *range.first, *range.second);
527 }
528 static_cast<VtkCompositeElementSelectionFilter*>(filter)
529 ->setSelectionArray("Selection", quality);
530 VtkVisPointSetItem* item =
531 new VtkVisPointSetItem(filter, parentItem, itemData);
532 this->addPipelineItem(item, this->createIndex(i, 0, item));
533 break;
534 }
535}
536
537void VtkVisPipeline::highlightGeoObject(const vtkPolyDataAlgorithm* source,
538 int index)
539{
541 int nSources = this->_rootItem->childCount();
542 for (int i = 0; i < nSources; i++)
543 {
544 auto* parentItem =
545 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
546 if (parentItem->algorithm() == source)
547 {
548 QList<QVariant> itemData;
549 itemData << "Selected GeoObject" << true;
550
551 VtkCompositeFilter* filter =
553 "VtkCompositeGeoObjectFilter",
554 parentItem->transformFilter());
555 static_cast<VtkCompositeGeoObjectFilter*>(filter)->SetIndex(index);
556 VtkVisPointSetItem* item =
557 new VtkVisPointSetItem(filter, parentItem, itemData);
558 QModelIndex parent_index =
559 static_cast<TreeModel*>(this)->index(i, 0, QModelIndex());
560 _highlighted_geo_index = this->addPipelineItem(item, parent_index);
561 break;
562 }
563 }
564}
565
567{
568 if (_highlighted_geo_index != QModelIndex())
569 {
571 _highlighted_geo_index = QModelIndex();
572 }
573}
574
576 vtkUnstructuredGridAlgorithm const* const source, unsigned index,
577 bool is_element)
578{
579 int nSources = this->_rootItem->childCount();
580 for (int i = 0; i < nSources; i++)
581 {
582 auto* parentItem =
583 static_cast<VtkVisPipelineItem*>(_rootItem->child(i));
584 if (parentItem->algorithm() == source)
585 {
586 QList<QVariant> itemData;
587 itemData << "Selected Mesh Component" << true;
588 QList<QVariant> selected_index;
589 selected_index << index << index;
590
591 VtkCompositeFilter* filter(nullptr);
592 if (is_element)
593 {
595 "VtkCompositeElementSelectionFilter",
596 parentItem->transformFilter());
597 static_cast<VtkCompositeElementSelectionFilter*>(filter)
598 ->setSelectionArray("vtkIdFilter_Ids");
599 static_cast<VtkCompositeElementSelectionFilter*>(filter)
600 ->SetUserVectorProperty("Threshold Between",
601 selected_index);
602 }
603 else
604 {
606 "VtkCompositeNodeSelectionFilter",
607 parentItem->transformFilter());
608 std::vector<unsigned> indices(1);
609 indices[0] = index;
610 static_cast<VtkCompositeNodeSelectionFilter*>(filter)
611 ->setSelectionArray(indices);
612 }
613 VtkVisPointSetItem* item =
614 new VtkVisPointSetItem(filter, parentItem, itemData);
615 QModelIndex parent_index =
616 static_cast<TreeModel*>(this)->index(i, 0, QModelIndex());
618 this->addPipelineItem(item, parent_index);
619 break;
620 }
621 }
622}
623
Definition of the GeoTreeModel class.
Definition of the LinearIntervalInterpolation class.
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
Definition of the MeshItem class.
Definition of the MeshModel class.
Definition of the Mesh class.
Definition of the StationTreeModel class.
Definition of the TreeModel class.
Definition of the VtkAlgorithmProperties class.
Definition of the VtkCompositeSelectionFilter class.
Definition of the VtkCompositeGeoObjectFilter class.
Definition of the VtkCompositeNodeSelectionFilter class.
Definition of the VtkFilterFactory class.
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
Definition of the VtkVisImageItem class.
Definition of the VtkVisPipelineItem class.
Definition of the VtkVisPipeline class.
Definition of the VtkVisPointSetItem class.
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:30
MeshLib::VtkMappedMeshSource * vtkSource() const
Returns the VTK object.
Definition MeshItem.h:44
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:28
virtual int childCount() const
Definition TreeItem.cpp:65
void appendChild(TreeItem *item)
Definition TreeItem.cpp:42
TreeItem * parentItem() const
Definition TreeItem.cpp:115
TreeItem * child(int row) const
Definition TreeItem.cpp:52
A hierarchical model for a tree implemented as a double-linked list.
Definition TreeModel.h:30
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:50
QModelIndex parent(const QModelIndex &index) const override
Definition TreeModel.cpp:81
TreeItem * _rootItem
Definition TreeModel.h:58
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:23
MeshQualityType
Describes a mesh quality metric.
Definition MeshEnums.h:70