Loading [MathJax]/extensions/tex2jax.js
OGS
VtkVisPipelineItem.cpp
Go to the documentation of this file.
1 
15 // ** INCLUDES **
16 #include "VtkVisPipelineItem.h"
17 
18 #include <vtkActor.h>
19 #include <vtkAlgorithm.h>
20 #include <vtkCellData.h>
21 #include <vtkDataSetMapper.h>
22 #include <vtkGenericDataObjectWriter.h>
23 #include <vtkImageActor.h>
24 #include <vtkPointData.h>
25 #include <vtkPointSet.h>
26 #include <vtkProperty.h>
27 #include <vtkRenderer.h>
28 #include <vtkSmartPointer.h>
29 #include <vtkTextureMapToPlane.h>
30 
31 #include <QMessageBox>
32 
33 #include "BaseLib/FileTools.h"
34 #include "QVtkDataSetMapper.h"
35 #include "VtkAlgorithmProperties.h"
36 #include "VtkCompositeFilter.h"
37 
39  vtkAlgorithm* algorithm, TreeItem* parentItem,
40  const QList<QVariant> data /*= QList<QVariant>()*/)
41  : TreeItem(data, parentItem),
42  _actor(nullptr),
43  _algorithm(algorithm),
44  _renderer(nullptr),
45  _compositeFilter(nullptr),
46  _vtkProps(nullptr)
47 {
48  auto* visParentItem = dynamic_cast<VtkVisPipelineItem*>(parentItem);
49  if (parentItem->parentItem())
50  {
51  _algorithm->SetInputConnection(
52  visParentItem->algorithm()->GetOutputPort());
53  }
54 }
55 
57  VtkCompositeFilter* compositeFilter, TreeItem* parentItem,
58  const QList<QVariant> data /*= QList<QVariant>()*/)
59  : TreeItem(data, parentItem),
60  _actor(nullptr),
61  _renderer(nullptr),
62  _compositeFilter(compositeFilter),
63  _vtkProps(nullptr)
64 {
66 }
67 
69 {
70  _renderer->RemoveActor(_actor);
71  _actor->Delete();
72  delete _compositeFilter;
73 }
74 
76 {
77  TreeItem* treeItem = TreeItem::child(row);
78  if (treeItem)
79  {
80  return dynamic_cast<VtkVisPipelineItem*>(treeItem);
81  }
82 
83  return nullptr;
84 }
85 
86 QVariant VtkVisPipelineItem::data(int column) const
87 {
88  if (column == 1)
89  {
90  return isVisible();
91  }
92 
93  return TreeItem::data(column);
94 }
95 
96 bool VtkVisPipelineItem::setData(int column, const QVariant& value)
97 {
98  if (column == 1)
99  {
100  setVisible(value.toBool());
101  return true;
102  }
103 
104  return TreeItem::setData(column, value);
105 }
107 {
108  return static_cast<bool>(_actor->GetVisibility());
109 }
110 
112 {
113  _actor->SetVisibility(static_cast<int>(visible));
114  _actor->Modified();
115  _renderer->Render();
116 }
117 
118 int VtkVisPipelineItem::writeToFile(const std::string& filename) const
119 {
120  if (!filename.empty())
121  {
122  return callVTKWriter(this->algorithm(), filename);
123  }
124  return 0;
125 }
126 
127 int VtkVisPipelineItem::callVTKWriter(vtkAlgorithm* algorithm,
128  const std::string& filename) const
129 {
130  // needs to be implemented in derived classes!
131  (void)algorithm;
132  (void)filename;
133  return 0;
134 }
135 
136 vtkProp3D* VtkVisPipelineItem::actor() const
137 {
138  return _actor;
139 }
140 
141 void VtkVisPipelineItem::setScale(double x, double y, double z) const
142 {
143  (void)x;
144  (void)y, (void)z;
145 }
146 
147 void VtkVisPipelineItem::setTranslation(double x, double y, double z) const
148 {
149  (void)x;
150  (void)y, (void)z;
151 }
152 
153 void VtkVisPipelineItem::setScaleOnChildren(double x, double y, double z) const
154 {
155  for (int i = 0; i < this->childCount(); ++i)
156  {
157  VtkVisPipelineItem* child = this->child(i);
158  child->setScale(x, y, z);
159  }
160 }
161 
163 {
164  // Reimplemented in subclass
165  (void)enable;
166 }
167 
169 {
170  for (int i = 0; i < this->childCount(); ++i)
171  {
172  VtkVisPipelineItem* child = this->child(i);
173  child->setBackfaceCulling(static_cast<int>(enable));
174  child->setBackfaceCullingOnChildren(static_cast<int>(enable));
175  }
176 }
177 
179 {
180  this->algorithm()->Update();
181  vtkDataSet* dataSet =
182  vtkDataSet::SafeDownCast(this->algorithm()->GetOutputDataObject(0));
183  QStringList dataSetAttributesList;
184  if (dataSet)
185  {
186  vtkPointData* pointData = dataSet->GetPointData();
187  // std::cout << " #point data arrays: " <<
188  // pointData->GetNumberOfArrays() << std::endl;
189  for (int i = 0; i < pointData->GetNumberOfArrays(); i++)
190  {
191  // std::cout << " Name: " << pointData->GetArrayName(i) <<
192  // std::endl;
193  dataSetAttributesList.push_back(QString("P-") +
194  pointData->GetArrayName(i));
195  }
196 
197  vtkCellData* cellData = dataSet->GetCellData();
198  // std::cout << " #cell data arrays: " << cellData->GetNumberOfArrays()
199  // << std::endl;
200  for (int i = 0; i < cellData->GetNumberOfArrays(); i++)
201  {
202  // std::cout << " Name: " << cellData->GetArrayName(i) <<
203  // std::endl;
204  dataSetAttributesList.push_back(QString("C-") +
205  cellData->GetArrayName(i));
206  }
207  }
208  return dataSetAttributesList;
209 }
Filename manipulation routines.
Definition of the QVtkDataSetMapper class.
Definition of the VtkAlgorithmProperties class.
Definition of the VtkCompositeFilter class.
Definition of the VtkVisPipelineItem class.
Objects nodes for the TreeModel.
Definition: TreeItem.h:28
virtual int childCount() const
Definition: TreeItem.cpp:65
TreeItem * parentItem() const
Definition: TreeItem.cpp:115
TreeItem * child(int row) const
Definition: TreeItem.cpp:52
int row() const
Definition: TreeItem.cpp:73
virtual bool setData(int column, const QVariant &value)
Definition: TreeItem.cpp:102
virtual QVariant data(int column) const
Definition: TreeItem.cpp:94
Is used to combine several filter in one VtkVisPipelineItem. You can use vtk filter and custom filter...
vtkAlgorithm * GetOutputAlgorithm() const
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...
virtual void setTranslation(double x, double y, double z) const
Translates the item in visualisation-space. This function is empty and needs to be implemented by der...
vtkProp3D * actor() const
Returns the actor as vtkProp3D.
QStringList getScalarArrayNames() const
Returns a list of array names prefixed with P- or C- for point and cell data.
vtkAlgorithm * algorithm() const
Returns the algorithm object.
virtual void setBackfaceCulling(bool enable) const
Enables / disables backface culling.
bool isVisible() const
Returns if the VTK object is visible in the visualization.
bool setData(int column, const QVariant &value) override
void setBackfaceCullingOnChildren(bool enable) const
Enables / disables backface culling on all children.
VtkVisPipelineItem(vtkAlgorithm *algorithm, TreeItem *parentItem, const QList< QVariant > data=QList< QVariant >())
Constructor for a source/filter object.
void setVisible(bool visible)
Sets the visibility of the VTK object in the visualization.
vtkAlgorithm * _algorithm
VtkCompositeFilter * _compositeFilter
VtkVisPipelineItem * child(int row) const
Returns a VtkVisPipelineItem.
virtual int callVTKWriter(vtkAlgorithm *algorithm, const std::string &filename) const
int writeToFile(const std::string &filename) const
Writes this algorithm's vtkDataSet (i.e. vtkPolyData or vtkUnstructuredGrid) to a vtk-file.
QVariant data(int column) const override
void setScaleOnChildren(double x, double y, double z) const
Sets the geometry and date scaling recursively on all children of this item.