OGS
VtkColorByHeightFilter.cpp
Go to the documentation of this file.
1
15// ** VTK INCLUDES **
17
18#include <vtkCellData.h>
19#include <vtkDoubleArray.h>
20#include <vtkInformation.h>
21#include <vtkInformationVector.h>
22#include <vtkLookupTable.h>
23#include <vtkObjectFactory.h>
24#include <vtkPointData.h>
25#include <vtkPolyData.h>
26#include <vtkSmartPointer.h>
27#include <vtkStreamingDemandDrivenPipeline.h>
28
29#include "BaseLib/Logging.h"
30#include "VtkColorLookupTable.h"
31
33
35 : ColorLookupTable(VtkColorLookupTable::New()), _tableRange{0.0, 0.0}
36{
37}
38
40
41void VtkColorByHeightFilter::PrintSelf(ostream& os, vtkIndent indent)
42{
43 this->Superclass::PrintSelf(os, indent);
44
45 double range[2];
46 ColorLookupTable->GetTableRange(range);
47 os << indent << "== VtkColorByHeightFilter ==" << endl;
48 os << indent << "Range: " << range[0] << "-" << range[1] << endl;
49 os << indent << "Interpolation Type:"
50 << static_cast<int>(ColorLookupTable->getInterpolationType()) << endl;
51}
52
54{
55 vtkMTimeType t1 = this->Superclass::GetMTime();
56 if (this->ColorLookupTable)
57 {
58 vtkMTimeType const t2 = this->ColorLookupTable->GetMTime();
59 return std::max(t1, t2);
60 }
61 return t1;
62}
63int VtkColorByHeightFilter::RequestData(vtkInformation* /*request*/,
64 vtkInformationVector** inputVector,
65 vtkInformationVector* outputVector)
66{
67 vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
68 vtkPolyData* input =
69 vtkPolyData::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT()));
70
71 vtkSmartPointer<vtkDoubleArray> colors =
72 vtkSmartPointer<vtkDoubleArray>::New();
73 colors->SetNumberOfComponents(1);
74 std::size_t nPoints = input->GetNumberOfPoints();
75 colors->SetNumberOfValues(nPoints);
76 colors->SetName("Colors");
77
78 // Inserts height values as a new scalar array
79 for (std::size_t i = 0; i < nPoints; i++)
80 {
81 double p[3];
82 input->GetPoint(i, p);
83 colors->SetValue(i, p[2]);
84 }
85
86 vtkInformation* outInfo = outputVector->GetInformationObject(0);
87 vtkPolyData* output =
88 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
89 output->CopyStructure(input);
90 output->GetPointData()->PassData(input->GetPointData());
91 output->GetCellData()->PassData(input->GetCellData());
92 output->GetPointData()->AddArray(colors);
93 output->GetPointData()->SetActiveScalars("Colors");
94
95 return 1;
96}
97
98void VtkColorByHeightFilter::SetTableRange(double min, double max)
99{
100 if (min < max)
101 {
102 this->_tableRange[0] = min;
103 this->_tableRange[1] = max;
104 this->ColorLookupTable->SetTableRange(min, max);
105 }
106 else
107 {
108 ERR("VtkColorByHeightFilter::SetLimits(min, max) - Limits not changed "
109 "because min value > max value.");
110 }
111}
112
114{
115 this->_tableRangeScaling = scale;
116 this->ColorLookupTable->SetTableRange(this->_tableRange[0] * scale,
117 this->_tableRange[1] * scale);
118}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
vtkStandardNewMacro(VtkColorByHeightFilter)
Definition of the VtkColorByHeightFilter class.
Definition of the VtkColorLookupTable class.
VTK filter object for colouring vtkPolyData objects based on z-coordinates.
void SetTableRangeScaling(double scale)
Sets the scaling of the color look-up table boundaries. This is used in VtkVisTabWidget when a parent...
~VtkColorByHeightFilter() override
vtkMTimeType GetMTime() override
This filter gets updated when the color look-up table was modified.
void PrintSelf(ostream &os, vtkIndent indent) override
Prints the mesh data to an output stream.
void SetTableRange(double min, double max)
Sets the boundaries for the color look-up table.
VtkColorLookupTable * ColorLookupTable
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
The filter logic.
Calculates and stores a colour lookup table.
DataHolderLib::LUTType getInterpolationType() const
Returns the type of interpolation used.