OGS
VtkCompositeThresholdFilter Class Reference

Detailed Description

Visualises only parts of meshes that are above/below/within given thresholds. In init() the threshold is first set to double min / max values. Set the threshold later on via SetUserVectorProperty() to the actual data range.

Definition at line 22 of file VtkCompositeThresholdFilter.h.

#include <VtkCompositeThresholdFilter.h>

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

Public Member Functions

 VtkCompositeThresholdFilter (vtkAlgorithm *inputAlgorithm)
 
 ~VtkCompositeThresholdFilter () override
 
void init () override
 
void SetUserProperty (QString name, QVariant value) override
 Sets a user property. This should be implemented by subclasses.
 
void SetUserVectorProperty (QString name, QList< QVariant > values) override
 Sets a vector user property. This should be implemented by subclasses.
 
- Public Member Functions inherited from VtkCompositeFilter
 VtkCompositeFilter (vtkAlgorithm *inputAlgorithm)
 Constructor.
 
 ~VtkCompositeFilter () override
 Destructor.
 
int GetInputDataObjectType () const
 
int GetOutputDataObjectType () const
 
vtkAlgorithm * GetOutputAlgorithm () const
 
- Public Member Functions inherited from VtkAlgorithmProperties
 VtkAlgorithmProperties (QObject *parent=nullptr)
 Constructor (sets default values)
 
 ~VtkAlgorithmProperties () override
 
vtkProperty * GetProperties () const
 Returns the vtk properties.
 
vtkTexture * GetTexture ()
 Returns a texture (if one has been assigned).
 
void SetTexture (vtkTexture *t)
 Sets a texture for the VtkVisPipelineItem.
 
vtkLookupTable * GetLookupTable (const QString &array_name)
 Returns the colour lookup table (if one has been assigned).
 
void RemoveLookupTable (const QString &array_name)
 Removes the lookup table for the given scalar.
 
void SetLookUpTable (const QString &array_name, vtkLookupTable *lut)
 Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
 
void SetLookUpTable (const QString &array_name, const QString &filename)
 Loads a predefined color lookup table from a file for the specified scalar array.
 
bool GetScalarVisibility () const
 Returns the scalar visibility.
 
void SetScalarVisibility (bool on)
 Sets the scalar visibility.
 
QString GetName () const
 Returns the name. This is set to the file path if it is a source algorithm.
 
void SetName (QString name)
 Sets the name.
 
bool IsRemovable () const
 Is this algorithm removable from the pipeline (view).
 
QMap< QString, QVariant > * GetAlgorithmUserProperties () const
 Returns a map of user properties.
 
QMap< QString, QList< QVariant > > * GetAlgorithmUserVectorProperties () const
 Returns a map of vector user properties.
 
QVariant GetUserProperty (QString name) const
 Returns the value of a user property.
 
QList< QVariant > GetUserVectorProperty (QString name) const
 Returns a list of values of a vector user property.
 
void SetActiveAttribute (QString name)
 Set the active attribute.
 
QString GetActiveAttribute () const
 Returns the desired active attribute.
 

Additional Inherited Members

- Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
 
- Protected Member Functions inherited from VtkCompositeFilter
double GetInitialRadius () const
 Calculates a 1/200th of the largest extension of the bounding box (this is used as default radius for various filters)
 
- Protected Attributes inherited from VtkCompositeFilter
int _inputDataObjectType
 
int _outputDataObjectType
 
vtkAlgorithm * _inputAlgorithm
 
vtkAlgorithm * _outputAlgorithm
 
- Protected Attributes inherited from VtkAlgorithmProperties
vtkProperty * _property
 
vtkTexture * _texture
 
bool _scalarVisibility
 
std::map< QString, vtkLookupTable * > _lut
 
QString _name
 
QString _activeAttributeName
 
bool _removable
 
QMap< QString, QVariant > * _algorithmUserProperties
 
QMap< QString, QList< QVariant > > * _algorithmUserVectorProperties
 

Constructor & Destructor Documentation

◆ VtkCompositeThresholdFilter()

VtkCompositeThresholdFilter::VtkCompositeThresholdFilter ( vtkAlgorithm * inputAlgorithm)
explicit

Definition at line 24 of file VtkCompositeThresholdFilter.cpp.

26 : VtkCompositeFilter(inputAlgorithm)
27{
28 this->init();
29}
VtkCompositeFilter(vtkAlgorithm *inputAlgorithm)
Constructor.

References init().

◆ ~VtkCompositeThresholdFilter()

VtkCompositeThresholdFilter::~VtkCompositeThresholdFilter ( )
overridedefault

Member Function Documentation

◆ init()

void VtkCompositeThresholdFilter::init ( )
overridevirtual

Implements VtkCompositeFilter.

Definition at line 33 of file VtkCompositeThresholdFilter.cpp.

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}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
vtkAlgorithm * _outputAlgorithm
vtkAlgorithm * _inputAlgorithm

References VtkCompositeFilter::_inputAlgorithm, VtkCompositeFilter::_inputDataObjectType, VtkCompositeFilter::_outputAlgorithm, VtkCompositeFilter::_outputDataObjectType, and WARN().

Referenced by VtkCompositeThresholdFilter().

◆ SetUserProperty()

void VtkCompositeThresholdFilter::SetUserProperty ( QString name,
QVariant value )
overridevirtual

Sets a user property. This should be implemented by subclasses.

Reimplemented from VtkAlgorithmProperties.

Definition at line 100 of file VtkCompositeThresholdFilter.cpp.

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}
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.

References VtkCompositeFilter::_outputAlgorithm, and VtkAlgorithmProperties::SetUserProperty().

◆ SetUserVectorProperty()

void VtkCompositeThresholdFilter::SetUserVectorProperty ( QString name,
QList< QVariant > values )
overridevirtual

Sets a vector user property. This should be implemented by subclasses.

Reimplemented from VtkAlgorithmProperties.

Definition at line 118 of file VtkCompositeThresholdFilter.cpp.

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}
virtual void SetUserVectorProperty(QString name, QList< QVariant > values)
Sets a vector user property. This should be implemented by subclasses.

References VtkCompositeFilter::_outputAlgorithm, and VtkAlgorithmProperties::SetUserVectorProperty().


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