OGS
VtkCompositeElementSelectionFilter.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
17
18#include <vtkDataSetSurfaceFilter.h>
19#include <vtkIdFilter.h>
20#include <vtkPointData.h>
21#include <vtkSmartPointer.h>
22#include <vtkThreshold.h>
23#include <vtkUnstructuredGrid.h>
24#include <vtkUnstructuredGridAlgorithm.h>
25
27#include "VtkColorLookupTable.h"
29#include "VtkPointsSource.h"
30
32 vtkAlgorithm* inputAlgorithm)
33 : VtkCompositeFilter(inputAlgorithm),
34 _range(0.0, 1.0),
35 _selection_name("Selection")
36{
37}
38
40{
41 double thresholdLower(_range.first);
42 double thresholdUpper(_range.second);
43 this->_inputDataObjectType = VTK_UNSTRUCTURED_GRID;
44 this->_outputDataObjectType = VTK_UNSTRUCTURED_GRID;
45
46 this->SetLookUpTable(QString::fromStdString(_selection_name),
47 this->GetLookupTable());
48 vtkSmartPointer<VtkAppendArrayFilter> selFilter(nullptr);
49 if (!_selection.empty())
50 {
51 selFilter = vtkSmartPointer<VtkAppendArrayFilter>::New();
52 selFilter->SetInputConnection(_inputAlgorithm->GetOutputPort());
53 selFilter->SetArray(_selection_name, _selection);
54 selFilter->Update();
55 }
56
57 vtkSmartPointer<vtkIdFilter> idFilter = vtkSmartPointer<vtkIdFilter>::New();
58 if (_selection.empty())
59 { // if the array is empty it is assumed that an existing array should be
60 // used
61 idFilter->SetInputConnection(_inputAlgorithm->GetOutputPort());
62 }
63 else
64 {
65 idFilter->SetInputConnection(selFilter->GetOutputPort());
66 }
67 idFilter->PointIdsOn();
68 idFilter->CellIdsOn();
69 idFilter->FieldDataOn();
70 idFilter->Update();
71
72 vtkThreshold* threshold = vtkThreshold::New();
73 threshold->SetInputConnection(idFilter->GetOutputPort());
74 threshold->SetInputArrayToProcess(0, 0, 0,
75 vtkDataObject::FIELD_ASSOCIATION_CELLS,
76 _selection_name.c_str());
77 threshold->SetSelectedComponent(0);
78 threshold->SetThresholdFunction(
79 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
80 threshold->SetLowerThreshold(thresholdLower);
81 threshold->SetUpperThreshold(thresholdUpper);
82 threshold->Update();
83
84 QList<QVariant> thresholdRangeList;
85 thresholdRangeList.push_back(thresholdLower);
86 thresholdRangeList.push_back(thresholdUpper);
87 (*_algorithmUserVectorProperties)["Threshold Between"] = thresholdRangeList;
88 _outputAlgorithm = threshold;
89}
90
92 const std::string& selection_name, const std::vector<double>& selection)
93{
94 _selection_name = selection_name;
95 _selection = selection;
96 init();
97}
98
100 QString name, QList<QVariant> values)
101{
103
104 if (name.compare("Threshold Between") == 0)
105 {
106 static_cast<vtkThreshold*>(_outputAlgorithm)
107 ->SetThresholdFunction(
108 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
109 static_cast<vtkThreshold*>(_outputAlgorithm)
110 ->SetLowerThreshold(values[0].toDouble());
111 static_cast<vtkThreshold*>(_outputAlgorithm)
112 ->SetUpperThreshold(values[1].toDouble());
113 }
114}
115
117{
119 lut->SetTableRange(0, 1);
120 DataHolderLib::Color a{{0, 0, 255, 255}}; // blue
121 DataHolderLib::Color b{{0, 255, 0, 255}}; // green
122 DataHolderLib::Color c{{255, 255, 0, 255}}; // yellow
123 DataHolderLib::Color d{{255, 0, 0, 255}}; // red
124 lut->setColor(1.0, a);
125 lut->setColor(0.5, b);
126 lut->setColor(0.25, c);
127 lut->setColor(0.1, d);
128 lut->Build();
129 return lut;
130}
Definition of the VtkAppendArrayFilter class.
Definition of the VtkColorLookupTable class.
Definition of the VtkCompositeSelectionFilter class.
Definition of the VtkCompositePointToGlyphFilter class.
Definition of the VtkPointsSource class.
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.
Calculates and stores a colour lookup table.
void setColor(double pos, DataHolderLib::Color const &color)
void Build() override
Builds the colour table based on the previously set parameters. This method should only be called aft...
static VtkColorLookupTable * New()
Create new objects with New() because of VTKs object reference counting.
VtkCompositeElementSelectionFilter(vtkAlgorithm *inputAlgorithm)
void setSelectionArray(const std::string &selection_name, const std::vector< double > &selection=std::vector< double >())
void SetUserVectorProperty(QString name, QList< QVariant > values) override
Sets a vector user property. This should be implemented by subclasses.
VtkColorLookupTable * GetLookupTable()
Returns a colour lookup table optimised for quality measures.
Is used to combine several filter in one VtkVisPipelineItem. You can use vtk filter and custom filter...
vtkAlgorithm * _outputAlgorithm
vtkAlgorithm * _inputAlgorithm
std::array< unsigned char, 4 > Color
Definition Color.h:24