OGS
VtkCompositeImageToCylindersFilter.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
17
18#include <vtkLookupTable.h>
19#include <vtkPointData.h>
20#include <vtkSmartPointer.h>
21#include <vtkTubeFilter.h>
22#include <vtkUnsignedCharArray.h>
23
24#include <QMap>
25#include <QString>
26#include <QVariant>
27#include <QVector>
28
30
32 vtkAlgorithm* inputAlgorithm)
33 : VtkCompositeFilter(inputAlgorithm), _lineFilter(nullptr)
34{
35 this->init();
36}
37
39{
40 this->_inputDataObjectType = VTK_IMAGE_DATA;
41 this->_outputDataObjectType = VTK_POLY_DATA;
42
44 _lineFilter->SetInputConnection(_inputAlgorithm->GetOutputPort());
45 _lineFilter->SetLengthScaleFactor(1);
46 (*_algorithmUserProperties)["LengthScaleFactor"] = 1.0;
47 _lineFilter->Update();
48
49 double range[2];
50 // The data is always on points
51 vtkDataSet::SafeDownCast(_lineFilter->GetOutputDataObject(0))
52 ->GetPointData()
53 ->GetScalars()
54 ->GetRange(range);
55
56 vtkLookupTable* colormap = vtkLookupTable::New();
57 colormap->SetTableRange(range[0], range[1]);
58 colormap->SetHueRange(0.0, 0.666);
59 colormap->SetNumberOfTableValues(256);
60 colormap->ForceBuild();
61 QList<QVariant> tableRangeList;
62 tableRangeList.push_back(range[0]);
63 tableRangeList.push_back(range[1]);
64 QList<QVariant> hueRangeList;
65 hueRangeList.push_back(0.0);
66 hueRangeList.push_back(0.666);
67 (*_algorithmUserVectorProperties)["TableRange"] = tableRangeList;
68 (*_algorithmUserVectorProperties)["HueRange"] = hueRangeList;
69
70 this->SetLookUpTable("P-Colors", colormap);
71
72 vtkTubeFilter* tubeFilter = vtkTubeFilter::New();
73 tubeFilter->SetInputConnection(_lineFilter->GetOutputPort());
74 tubeFilter->CappingOn();
75 tubeFilter->SetNumberOfSides(6);
76 tubeFilter->SetRadius(_lineFilter->GetImageSpacing() * 0.25);
77 (*_algorithmUserProperties)["NumberOfColors"] = 256;
78 (*_algorithmUserProperties)["Capping"] = true;
79 (*_algorithmUserProperties)["NumberOfSides"] = 6;
80 (*_algorithmUserProperties)["RadiusFactor"] = 0.25;
81
82 _outputAlgorithm = tubeFilter;
83}
84
86 QVariant value)
87{
89
90 _lineFilter->SetUserProperty(name, value);
91
92 // VtkImageDataToLinePolyDataFilter is equal to _firstAlgorithm
93 // vtkTubeFilter is equal _outputAlgorithm
94 if (name.compare("NumberOfColors") == 0)
95 {
96 vtkLookupTable* lut = this->GetLookupTable("P-Colors");
97 if (lut)
98 {
99 lut->SetNumberOfTableValues(value.toInt());
100 }
101 }
102 else if (name.compare("NumberOfSides") == 0)
103 {
104 static_cast<vtkTubeFilter*>(_outputAlgorithm)
105 ->SetNumberOfSides(value.toInt());
106 }
107 else if (name.compare("Capping") == 0)
108 {
109 static_cast<vtkTubeFilter*>(_outputAlgorithm)
110 ->SetCapping(value.toBool());
111 }
112 else if (name.compare("RadiusFactor") == 0)
113 {
114 static_cast<vtkTubeFilter*>(_outputAlgorithm)
115 ->SetRadius(_lineFilter->GetImageSpacing() * value.toDouble());
116 }
117}
118
120 QString name, QList<QVariant> values)
121{
123
124 _lineFilter->SetUserVectorProperty(name, values);
125
126 if (name.compare("TableRange") == 0)
127 {
128 vtkLookupTable* lut = this->GetLookupTable("P-Colors");
129 if (lut)
130 {
131 lut->SetTableRange(values[0].toDouble(), values[1].toDouble());
132 }
133 }
134 else if (name.compare("HueRange") == 0)
135 {
136 vtkLookupTable* lut = this->GetLookupTable("P-Colors");
137 if (lut)
138 {
139 lut->SetHueRange(values[0].toDouble(), values[1].toDouble());
140 }
141 }
142}
143
Definition of the VtkCompositeImageToCylindersFilter class.
Definition of the VtkImageDataToLinePolyDataFilter class.
vtkLookupTable * GetLookupTable(const QString &array_name)
Returns the colour lookup table (if one has been assigned).
virtual void SetUserVectorProperty(QString name, QList< QVariant > values)
Sets a vector user property. This should be implemented by subclasses.
void SetLookUpTable(const QString &array_name, vtkLookupTable *lut)
Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.
Is used to combine several filter in one VtkVisPipelineItem. You can use vtk filter and custom filter...
vtkAlgorithm * _outputAlgorithm
vtkAlgorithm * _inputAlgorithm
void SetUserVectorProperty(QString name, QList< QVariant > values) override
Sets a vector user property. This should be implemented by subclasses.
void SetUserProperty(QString name, QVariant value) override
Sets a user property. This should be implemented by subclasses.
VtkCompositeImageToCylindersFilter(vtkAlgorithm *inputAlgorithm)
static VtkImageDataToLinePolyDataFilter * New()
Create new objects with New() because of VTKs reference counting.
void SetUserProperty(QString name, QVariant value) override
Sets a user property.