OGS
VtkCompositeElementSelectionFilter.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4// ** INCLUDES **
6
7#include <vtkDataSetSurfaceFilter.h>
8#include <vtkIdFilter.h>
9#include <vtkPointData.h>
10#include <vtkSmartPointer.h>
11#include <vtkThreshold.h>
12#include <vtkUnstructuredGrid.h>
13#include <vtkUnstructuredGridAlgorithm.h>
14
16#include "VtkColorLookupTable.h"
18#include "VtkPointsSource.h"
19
21 vtkAlgorithm* inputAlgorithm)
22 : VtkCompositeFilter(inputAlgorithm),
23 _range(0.0, 1.0),
24 _selection_name("Selection")
25{
26}
27
29{
30 double thresholdLower(_range.first);
31 double thresholdUpper(_range.second);
32 this->_inputDataObjectType = VTK_UNSTRUCTURED_GRID;
33 this->_outputDataObjectType = VTK_UNSTRUCTURED_GRID;
34
35 this->SetLookUpTable(QString::fromStdString(_selection_name),
36 this->GetLookupTable());
37 vtkSmartPointer<VtkAppendArrayFilter> selFilter(nullptr);
38 if (!_selection.empty())
39 {
40 selFilter = vtkSmartPointer<VtkAppendArrayFilter>::New();
41 selFilter->SetInputConnection(_inputAlgorithm->GetOutputPort());
42 selFilter->SetArray(_selection_name, _selection);
43 selFilter->Update();
44 }
45
46 vtkSmartPointer<vtkIdFilter> idFilter = vtkSmartPointer<vtkIdFilter>::New();
47 if (_selection.empty())
48 { // if the array is empty it is assumed that an existing array should be
49 // used
50 idFilter->SetInputConnection(_inputAlgorithm->GetOutputPort());
51 }
52 else
53 {
54 idFilter->SetInputConnection(selFilter->GetOutputPort());
55 }
56 idFilter->PointIdsOn();
57 idFilter->CellIdsOn();
58 idFilter->FieldDataOn();
59 idFilter->Update();
60
61 vtkThreshold* threshold = vtkThreshold::New();
62 threshold->SetInputConnection(idFilter->GetOutputPort());
63 threshold->SetInputArrayToProcess(0, 0, 0,
64 vtkDataObject::FIELD_ASSOCIATION_CELLS,
65 _selection_name.c_str());
66 threshold->SetSelectedComponent(0);
67 threshold->SetThresholdFunction(
68 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
69 threshold->SetLowerThreshold(thresholdLower);
70 threshold->SetUpperThreshold(thresholdUpper);
71 threshold->Update();
72
73 QList<QVariant> thresholdRangeList;
74 thresholdRangeList.push_back(thresholdLower);
75 thresholdRangeList.push_back(thresholdUpper);
76 (*_algorithmUserVectorProperties)["Threshold Between"] = thresholdRangeList;
77 _outputAlgorithm = threshold;
78}
79
81 const std::string& selection_name, const std::vector<double>& selection)
82{
83 _selection_name = selection_name;
84 _selection = selection;
85 init();
86}
87
89 QString name, QList<QVariant> values)
90{
92
93 if (name.compare("Threshold Between") == 0)
94 {
95 static_cast<vtkThreshold*>(_outputAlgorithm)
96 ->SetThresholdFunction(
97 vtkThreshold::ThresholdType::THRESHOLD_BETWEEN);
98 static_cast<vtkThreshold*>(_outputAlgorithm)
99 ->SetLowerThreshold(values[0].toDouble());
100 static_cast<vtkThreshold*>(_outputAlgorithm)
101 ->SetUpperThreshold(values[1].toDouble());
102 }
103}
104
106{
108 lut->SetTableRange(0, 1);
109 DataHolderLib::Color a{{0, 0, 255, 255}}; // blue
110 DataHolderLib::Color b{{0, 255, 0, 255}}; // green
111 DataHolderLib::Color c{{255, 255, 0, 255}}; // yellow
112 DataHolderLib::Color d{{255, 0, 0, 255}}; // red
113 lut->setColor(1.0, a);
114 lut->setColor(0.5, b);
115 lut->setColor(0.25, c);
116 lut->setColor(0.1, d);
117 lut->Build();
118 return lut;
119}
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.
vtkAlgorithm * _outputAlgorithm
VtkCompositeFilter(vtkAlgorithm *inputAlgorithm)
Constructor.
vtkAlgorithm * _inputAlgorithm
std::array< unsigned char, 4 > Color
Definition Color.h:12