OGS
VtkCompositeThresholdFilter.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
17
18#include <vtkCellData.h>
19#include <vtkThreshold.h>
20#include <vtkUnstructuredGrid.h>
21
22#include "BaseLib/Logging.h"
23
25 vtkAlgorithm* inputAlgorithm)
26 : VtkCompositeFilter(inputAlgorithm)
27{
28 this->init();
29}
30
32
34{
35 // Set meta information about input and output data types
36 this->_inputDataObjectType = VTK_DATA_SET;
37 this->_outputDataObjectType = VTK_UNSTRUCTURED_GRID;
38
39 // Because this is the only filter here we cannot use vtkSmartPointer
40 vtkThreshold* threshold = vtkThreshold::New();
41 threshold->SetInputConnection(_inputAlgorithm->GetOutputPort());
42
43 // Use first array of parent as input array
44 _inputAlgorithm->Update();
45 vtkDataSet* dataSet =
46 vtkDataSet::SafeDownCast(_inputAlgorithm->GetOutputDataObject(0));
47 vtkDataSetAttributes* pointAttributes =
48 dataSet->GetAttributes(vtkDataObject::AttributeTypes::POINT);
49 vtkDataSetAttributes* cellAttributes =
50 dataSet->GetAttributes(vtkDataObject::AttributeTypes::CELL);
51 if (pointAttributes->GetNumberOfArrays() > 0)
52 {
53 threshold->SetInputArrayToProcess(
54 0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS,
55 pointAttributes->GetArray(0)->GetName());
56 }
57 else if (cellAttributes->GetNumberOfArrays() > 0)
58 {
59 threshold->SetInputArrayToProcess(
60 0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS,
61 cellAttributes->GetArray(0)->GetName());
62 }
63 else
64 {
65 WARN("Threshold filter could not find an array on its input object!");
66 return;
67 }
68
69 // Sets a filter property which will be user editable
70 threshold->SetSelectedComponent(0);
71
72 // Setting the threshold to min / max values to ensure that the whole data
73 // is first processed. This is needed for correct lookup table generation.
74 const double dMin = std::numeric_limits<double>::lowest();
75 const double dMax = std::numeric_limits<double>::max();
76 threshold->SetThresholdFunction(
77 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
78 threshold->SetLowerThreshold(dMin);
79 threshold->SetUpperThreshold(dMax);
80
81 // Create a list for the ThresholdBetween (vector) property.
82 QList<QVariant> thresholdRangeList;
83 // Insert the values (same values as above)
84 thresholdRangeList.push_back(dMin);
85 thresholdRangeList.push_back(dMax);
86 // Put that list in the property map
87 (*_algorithmUserVectorProperties)["Range"] = thresholdRangeList;
88
89 // Make a new entry in the property map for the SelectedComponent property
90 (*_algorithmUserProperties)["Selected Component"] = 0;
91
92 // Must all scalars match the criterium
93 threshold->SetAllScalars(1);
94 (*_algorithmUserProperties)["Evaluate all points"] = true;
95
96 // The threshold filter is last one and so it is also the _outputAlgorithm
97 _outputAlgorithm = threshold;
98}
99
100void VtkCompositeThresholdFilter::SetUserProperty(QString name, QVariant value)
101{
103
104 // Use the same name as in init()
105 if (name.compare("Selected Component") == 0)
106 {
107 // Set the property on the algorithm
108 static_cast<vtkThreshold*>(_outputAlgorithm)
109 ->SetSelectedComponent(value.toInt());
110 }
111 else if (name.compare("Evaluate all points") == 0)
112 {
113 static_cast<vtkThreshold*>(_outputAlgorithm)
114 ->SetAllScalars(value.toBool());
115 }
116}
117
119 QList<QVariant> values)
120{
122
123 // Use the same name as in init()
124 if (name.compare("Range") == 0)
125 {
126 // Set the vector property on the algorithm
127 static_cast<vtkThreshold*>(_outputAlgorithm)
128 ->SetThresholdFunction(
129 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
130 static_cast<vtkThreshold*>(_outputAlgorithm)
131 ->SetLowerThreshold(values[0].toDouble());
132 static_cast<vtkThreshold*>(_outputAlgorithm)
133 ->SetUpperThreshold(values[1].toDouble());
134 }
135}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Definition of the VtkCompositeThresholdFilter class.
virtual void SetUserVectorProperty(QString name, QList< QVariant > values)
Sets a vector user property. This should be implemented by subclasses.
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
VtkCompositeThresholdFilter(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.
~VtkCompositeThresholdFilter() override