OGS
VtkVisPipeline Class Reference

Detailed Description

VtkVisPipeline manages the VTK visualization. It is a TreeModel and provides functions for adding and removing OGS Model or vtkAlgorithm objects.

Definition at line 44 of file VtkVisPipeline.h.

#include <VtkVisPipeline.h>

Inheritance diagram for VtkVisPipeline:
[legend]
Collaboration diagram for VtkVisPipeline:
[legend]

Public Slots

void addPipelineItem (MeshModel *model, const QModelIndex &idx)
 Adds the given Model to the pipeline.
void addPipelineItem (GeoTreeModel *model, const std::string &name, GeoLib::GEOTYPE type)
void addPipelineItem (StationTreeModel *model, const std::string &name)
QModelIndex addPipelineItem (VtkVisPipelineItem *item, const QModelIndex &parent)
QModelIndex addPipelineItem (vtkAlgorithm *source, QModelIndex parent=QModelIndex())
 Inserts the vtkAlgorithm as a child of the given QModelIndex to the pipeline.
void removeSourceItem (MeshModel *model, const QModelIndex &idx)
 Removes the given Model (and all attached vtkAlgorithms) from the pipeline.
void removeSourceItem (GeoTreeModel *model, const std::string &name, GeoLib::GEOTYPE type)
void removeSourceItem (StationTreeModel *model, const std::string &name)
void removePipelineItem (QModelIndex index)
 Removes the vtkAlgorithm at the given QModelIndex (and all attached vtkAlgorithms) from the pipeline.
void highlightGeoObject (const vtkPolyDataAlgorithm *source, int index)
 Applies a VtkCompositeGeoObjectFilter to add a specific index of the given geometry-source to the pipeline for highlighted display in the render window.
void removeHighlightedGeoObject ()
 Removes the currently highlighted geometry-object.
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 pipeline for highlighted display in the render window.
void removeHighlightedMeshComponent ()
 Removes the currently highlighted mesh component.
Public Slots inherited from TreeModel
void updateData ()

Signals

void vtkVisPipelineChanged () const
 Is emitted when a pipeline item was added or removed.
void itemSelected (const QModelIndex &) const

Public Member Functions

 VtkVisPipeline (vtkRenderer *renderer, QObject *parent=nullptr)
bool setData (const QModelIndex &index, const QVariant &value, int role=Qt::EditRole) override
 Emits vtkVisPipelineChanged() and calls base class method.
void addLight (const GeoLib::Point &pos)
 Adds a light to the scene at the given coordinates.
vtkLight * getLight (const GeoLib::Point &pos) const
 Returns a light (or nullptr) for the given coordinates.
void removeLight (const GeoLib::Point &pos)
 Removes a light at the given coordinates (if possible).
QColor getBGColor () const
 Returns the background-colour of the scene.
void setBGColor (const QColor &color)
 Sets the background-colour of the scene.
QModelIndex getIndex (vtkProp3D *actor)
 Returns the QModelIndex of VtkVisPipelineItem which actor is the given one.
Qt::ItemFlags flags (const QModelIndex &index) const override
void loadFromFile (QString filename)
 Loads a vtk object from the given file and adds it to the pipeline.
void resetCameraOnAddOrRemove (bool reset)
 Defaults to on.
void setGlobalSuperelevation (double factor) const
 Sets a global superelevation factor on all source items and resets the factor on other items to 1.
void setGlobalBackfaceCulling (bool enable) const
 Enables / disables backface culling on all actors.
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.
Public Member Functions inherited from TreeModel
 TreeModel (QObject *parent=nullptr)
 ~TreeModel () override
QVariant data (const QModelIndex &index, int role) const override
bool setData (const QModelIndex &index, const QVariant &value, int role) override
Qt::ItemFlags flags (const QModelIndex &index) const override
TreeItemgetItem (const QModelIndex &index) const
QVariant headerData (int section, Qt::Orientation orientation, int role=Qt::DisplayRole) const override
QModelIndex index (int row, int column, const QModelIndex &parent=QModelIndex()) const override
QModelIndex parent (const QModelIndex &index) const override
bool removeRows (int position, int count, const QModelIndex &parent) override
int rowCount (const QModelIndex &parent=QModelIndex()) const override
int columnCount (const QModelIndex &parent=QModelIndex()) const override
TreeItemrootItem () const

Private Member Functions

void listArrays (vtkDataSet *dataSet)

Private Attributes

vtkRenderer * _renderer
QVector< vtkAlgorithm * > _sources
std::list< vtkLight * > _lights
QMap< vtkProp3D *, QModelIndex > _actorMap
bool _resetCameraOnAddOrRemove
QModelIndex _highlighted_geo_index
QModelIndex _highlighted_mesh_component

Additional Inherited Members

Protected Attributes inherited from TreeModel
TreeItem_rootItem

Constructor & Destructor Documentation

◆ VtkVisPipeline()

VtkVisPipeline::VtkVisPipeline ( vtkRenderer * renderer,
QObject * parent = nullptr )
explicit

Definition at line 53 of file VtkVisPipeline.cpp.

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}
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
void setBGColor(const QColor &color)
Sets the background-colour of the scene.
bool _resetCameraOnAddOrRemove
vtkRenderer * _renderer

References TreeModel::TreeModel(), _renderer, _resetCameraOnAddOrRemove, TreeModel::_rootItem, TreeModel::parent(), and setBGColor().

Member Function Documentation

◆ addLight()

void VtkVisPipeline::addLight ( const GeoLib::Point & pos)

Adds a light to the scene at the given coordinates.

Definition at line 81 of file VtkVisPipeline.cpp.

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}
std::list< vtkLight * > _lights

References _lights, and _renderer.

◆ addPipelineItem [1/5]

void VtkVisPipeline::addPipelineItem ( GeoTreeModel * model,
const std::string & name,
GeoLib::GEOTYPE type )
slot

Definition at line 275 of file VtkVisPipeline.cpp.

278{
279 addPipelineItem(model->vtkSource(name, type));
280}
vtkPolyDataAlgorithm * vtkSource(const std::string &name, GeoLib::GEOTYPE type) const
Returns the vtk-object indicated by type of the geometry indicated by name.
void addPipelineItem(MeshModel *model, const QModelIndex &idx)
Adds the given Model to the pipeline.

References addPipelineItem(), and GeoTreeModel::vtkSource().

◆ addPipelineItem [2/5]

void VtkVisPipeline::addPipelineItem ( MeshModel * model,
const QModelIndex & idx )
slot

Adds the given Model to the pipeline.

Definition at line 288 of file VtkVisPipeline.cpp.

289{
290 addPipelineItem(static_cast<MeshItem*>(model->getItem(idx))->vtkSource());
291}
TreeItem * getItem(const QModelIndex &index) const

References addPipelineItem(), TreeModel::getItem(), and MeshItem::vtkSource().

Referenced by addPipelineItem(), addPipelineItem(), addPipelineItem(), addPipelineItem(), highlightGeoObject(), highlightMeshComponent(), loadFromFile(), and showMeshElementQuality().

◆ addPipelineItem [3/5]

void VtkVisPipeline::addPipelineItem ( StationTreeModel * model,
const std::string & name )
slot

Definition at line 282 of file VtkVisPipeline.cpp.

284{
285 addPipelineItem(model->vtkSource(name));
286}
vtkPolyDataAlgorithm * vtkSource(const std::string &name) const

References addPipelineItem(), and StationTreeModel::vtkSource().

◆ addPipelineItem [4/5]

QModelIndex VtkVisPipeline::addPipelineItem ( vtkAlgorithm * source,
QModelIndex parent = QModelIndex() )
slot

Inserts the vtkAlgorithm as a child of the given QModelIndex to the pipeline.

Definition at line 335 of file VtkVisPipeline.cpp.

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}

References addPipelineItem(), TreeModel::getItem(), and TreeModel::parent().

◆ addPipelineItem [5/5]

QModelIndex VtkVisPipeline::addPipelineItem ( VtkVisPipelineItem * item,
const QModelIndex & parent )
slot

Definition at line 293 of file VtkVisPipeline.cpp.

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}
virtual int childCount() const
Definition TreeItem.cpp:54
void appendChild(TreeItem *item)
Definition TreeItem.cpp:31
TreeItem * parentItem() const
Definition TreeItem.cpp:104
QModelIndex index(int row, int column, const QModelIndex &parent=QModelIndex()) const override
Definition TreeModel.cpp:39
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 Initialize(vtkRenderer *renderer)=0
Initializes vtkMapper and vtkActor necessary for visualization of the item and sets the item's proper...
void itemSelected(const QModelIndex &) const
void vtkVisPipelineChanged() const
Is emitted when a pipeline item was added or removed.
QMap< vtkProp3D *, QModelIndex > _actorMap

References _actorMap, _renderer, _resetCameraOnAddOrRemove, TreeModel::_rootItem, VtkVisPipelineItem::actor(), VtkVisPipelineItem::algorithm(), TreeItem::appendChild(), TreeItem::childCount(), TreeModel::index(), VtkVisPipelineItem::Initialize(), itemSelected(), TreeModel::parent(), TreeItem::parentItem(), VtkVisPipelineItem::setScale(), and vtkVisPipelineChanged().

◆ flags()

Qt::ItemFlags VtkVisPipeline::flags ( const QModelIndex & index) const
override

Definition at line 152 of file VtkVisPipeline.cpp.

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}

References TreeModel::index().

◆ getBGColor()

QColor VtkVisPipeline::getBGColor ( ) const

Returns the background-colour of the scene.

Definition at line 131 of file VtkVisPipeline.cpp.

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}

References _renderer.

◆ getIndex()

QModelIndex VtkVisPipeline::getIndex ( vtkProp3D * actor)

Returns the QModelIndex of VtkVisPipelineItem which actor is the given one.

Definition at line 147 of file VtkVisPipeline.cpp.

148{
149 return _actorMap.value(actor, QModelIndex());
150}

References _actorMap.

◆ getLight()

vtkLight * VtkVisPipeline::getLight ( const GeoLib::Point & pos) const

Returns a light (or nullptr) for the given coordinates.

Definition at line 99 of file VtkVisPipeline.cpp.

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}

References _lights.

◆ highlightGeoObject

void VtkVisPipeline::highlightGeoObject ( const vtkPolyDataAlgorithm * source,
int index )
slot

Applies a VtkCompositeGeoObjectFilter to add a specific index of the given geometry-source to the pipeline for highlighted display in the render window.

Definition at line 526 of file VtkVisPipeline.cpp.

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}
static VtkCompositeFilter * CreateCompositeFilter(QString type, vtkAlgorithm *inputAlgorithm)
Creates a composite filter by name.
void removeHighlightedGeoObject()
Removes the currently highlighted geometry-object.
QModelIndex _highlighted_geo_index
constexpr List filter(Pred, SomeListOfTypes<> *)
Definition TMP.h:42

References TreeModel::TreeModel(), _highlighted_geo_index, TreeModel::_rootItem, addPipelineItem(), VtkFilterFactory::CreateCompositeFilter(), TreeModel::index(), and removeHighlightedGeoObject().

◆ highlightMeshComponent

void VtkVisPipeline::highlightMeshComponent ( vtkUnstructuredGridAlgorithm const *const source,
unsigned index,
bool is_element )
slot

Applies a VtkCompositeSelectionFilter to add a specific component of the given mesh-source to the pipeline for highlighted display in the render window.

Definition at line 564 of file VtkVisPipeline.cpp.

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}
QModelIndex _highlighted_mesh_component

References TreeModel::TreeModel(), _highlighted_mesh_component, TreeModel::_rootItem, addPipelineItem(), VtkFilterFactory::CreateCompositeFilter(), and TreeModel::index().

◆ itemSelected

void VtkVisPipeline::itemSelected ( const QModelIndex & ) const
signal

Referenced by addPipelineItem().

◆ listArrays()

void VtkVisPipeline::listArrays ( vtkDataSet * dataSet)
private

Definition at line 465 of file VtkVisPipeline.cpp.

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}
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

References ERR(), and INFO().

Referenced by loadFromFile().

◆ loadFromFile()

void VtkVisPipeline::loadFromFile ( QString filename)

Loads a vtk object from the given file and adds it to the pipeline.

Definition at line 167 of file VtkVisPipeline.cpp.

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}
void listArrays(vtkDataSet *dataSet)

References addPipelineItem(), ERR(), INFO(), and listArrays().

◆ removeHighlightedGeoObject

void VtkVisPipeline::removeHighlightedGeoObject ( )
slot

Removes the currently highlighted geometry-object.

Definition at line 555 of file VtkVisPipeline.cpp.

556{
557 if (_highlighted_geo_index != QModelIndex())
558 {
560 _highlighted_geo_index = QModelIndex();
561 }
562}
void removePipelineItem(QModelIndex index)
Removes the vtkAlgorithm at the given QModelIndex (and all attached vtkAlgorithms) from the pipeline.

References _highlighted_geo_index, and removePipelineItem().

Referenced by highlightGeoObject().

◆ removeHighlightedMeshComponent

void VtkVisPipeline::removeHighlightedMeshComponent ( )
slot

Removes the currently highlighted mesh component.

Definition at line 613 of file VtkVisPipeline.cpp.

614{
615 if (_highlighted_mesh_component != QModelIndex())
616 {
618 _highlighted_mesh_component = QModelIndex();
619 }
620}

References _highlighted_mesh_component, and removePipelineItem().

◆ removeLight()

void VtkVisPipeline::removeLight ( const GeoLib::Point & pos)

Removes a light at the given coordinates (if possible).

Definition at line 114 of file VtkVisPipeline.cpp.

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}

References _lights, and _renderer.

◆ removePipelineItem

void VtkVisPipeline::removePipelineItem ( QModelIndex index)
slot

Removes the vtkAlgorithm at the given QModelIndex (and all attached vtkAlgorithms) from the pipeline.

Definition at line 436 of file VtkVisPipeline.cpp.

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}
bool removeRows(int position, int count, const QModelIndex &parent) override

References _actorMap, _renderer, _resetCameraOnAddOrRemove, TreeModel::index(), TreeModel::removeRows(), and vtkVisPipelineChanged().

Referenced by removeHighlightedGeoObject(), removeHighlightedMeshComponent(), removeSourceItem(), removeSourceItem(), and removeSourceItem().

◆ removeSourceItem [1/3]

void VtkVisPipeline::removeSourceItem ( GeoTreeModel * model,
const std::string & name,
GeoLib::GEOTYPE type )
slot

Definition at line 389 of file VtkVisPipeline.cpp.

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}

References TreeModel::_rootItem, VtkVisPipelineItem::algorithm(), TreeModel::getItem(), TreeModel::index(), removePipelineItem(), and GeoTreeModel::vtkSource().

◆ removeSourceItem [2/3]

void VtkVisPipeline::removeSourceItem ( MeshModel * model,
const QModelIndex & idx )
slot

Removes the given Model (and all attached vtkAlgorithms) from the pipeline.

Definition at line 420 of file VtkVisPipeline.cpp.

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}

References TreeModel::_rootItem, VtkVisPipelineItem::algorithm(), TreeModel::getItem(), TreeModel::index(), and removePipelineItem().

◆ removeSourceItem [3/3]

void VtkVisPipeline::removeSourceItem ( StationTreeModel * model,
const std::string & name )
slot

Definition at line 405 of file VtkVisPipeline.cpp.

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}

References TreeModel::_rootItem, VtkVisPipelineItem::algorithm(), TreeModel::getItem(), TreeModel::index(), removePipelineItem(), and StationTreeModel::vtkSource().

◆ resetCameraOnAddOrRemove()

void VtkVisPipeline::resetCameraOnAddOrRemove ( bool reset)
inline

Defaults to on.

Definition at line 80 of file VtkVisPipeline.h.

References _resetCameraOnAddOrRemove.

◆ setBGColor()

void VtkVisPipeline::setBGColor ( const QColor & color)

Sets the background-colour of the scene.

Definition at line 140 of file VtkVisPipeline.cpp.

141{
142 QSettings settings;
143 settings.setValue("VtkBackgroundColor", color);
144 _renderer->SetBackground(color.redF(), color.greenF(), color.blueF());
145}

References _renderer.

Referenced by VtkVisPipeline().

◆ setData()

bool VtkVisPipeline::setData ( const QModelIndex & index,
const QVariant & value,
int role = Qt::EditRole )
override

Emits vtkVisPipelineChanged() and calls base class method.

Definition at line 73 of file VtkVisPipeline.cpp.

75{
77
78 return TreeModel::setData(index, value, role);
79}
bool setData(const QModelIndex &index, const QVariant &value, int role) override

References TreeModel::index(), TreeModel::setData(), and vtkVisPipelineChanged().

◆ setGlobalBackfaceCulling()

void VtkVisPipeline::setGlobalBackfaceCulling ( bool enable) const

Enables / disables backface culling on all actors.

Definition at line 260 of file VtkVisPipeline.cpp.

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}

References TreeModel::_rootItem, VtkVisPipelineItem::setBackfaceCulling(), and vtkVisPipelineChanged().

◆ setGlobalSuperelevation()

void VtkVisPipeline::setGlobalSuperelevation ( double factor) const

Sets a global superelevation factor on all source items and resets the factor on other items to 1.

Definition at line 245 of file VtkVisPipeline.cpp.

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}

References TreeModel::_rootItem, VtkVisPipelineItem::setScale(), and vtkVisPipelineChanged().

◆ showMeshElementQuality()

void VtkVisPipeline::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.

Definition at line 483 of file VtkVisPipeline.cpp.

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}
std::string MeshQualityType2String(const MeshQualityType t)

References TreeModel::_rootItem, addPipelineItem(), VtkFilterFactory::CreateCompositeFilter(), and MeshLib::ELEMENTSIZE.

◆ vtkVisPipelineChanged

void VtkVisPipeline::vtkVisPipelineChanged ( ) const
signal

Is emitted when a pipeline item was added or removed.

Referenced by addPipelineItem(), removePipelineItem(), setData(), setGlobalBackfaceCulling(), and setGlobalSuperelevation().

Member Data Documentation

◆ _actorMap

QMap<vtkProp3D*, QModelIndex> VtkVisPipeline::_actorMap
private

Definition at line 132 of file VtkVisPipeline.h.

Referenced by addPipelineItem(), getIndex(), and removePipelineItem().

◆ _highlighted_geo_index

QModelIndex VtkVisPipeline::_highlighted_geo_index
private

Definition at line 135 of file VtkVisPipeline.h.

Referenced by highlightGeoObject(), and removeHighlightedGeoObject().

◆ _highlighted_mesh_component

QModelIndex VtkVisPipeline::_highlighted_mesh_component
private

Definition at line 136 of file VtkVisPipeline.h.

Referenced by highlightMeshComponent(), and removeHighlightedMeshComponent().

◆ _lights

std::list<vtkLight*> VtkVisPipeline::_lights
private

Definition at line 131 of file VtkVisPipeline.h.

Referenced by addLight(), getLight(), and removeLight().

◆ _renderer

vtkRenderer* VtkVisPipeline::_renderer
private

◆ _resetCameraOnAddOrRemove

bool VtkVisPipeline::_resetCameraOnAddOrRemove
private

◆ _sources

QVector<vtkAlgorithm*> VtkVisPipeline::_sources
private

Definition at line 130 of file VtkVisPipeline.h.


The documentation for this class was generated from the following files: