OGS
VtkTextureOnSurfaceFilter.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
17
18#include <vtkCellData.h>
19#include <vtkFloatArray.h>
20#include <vtkImageAlgorithm.h>
21#include <vtkImageShiftScale.h>
22#include <vtkInformation.h>
23#include <vtkInformationVector.h>
24#include <vtkObjectFactory.h>
25#include <vtkPointData.h>
26#include <vtkProperty.h>
27#include <vtkSmartPointer.h>
28#include <vtkStreamingDemandDrivenPipeline.h>
29#include <vtkTexture.h>
30
31#include "BaseLib/Logging.h"
32#include "MathLib/MathTools.h"
33#include "VtkVisHelper.h"
34
36
39
40void VtkTextureOnSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
41{
42 this->Superclass::PrintSelf(os, indent);
43}
44
45int VtkTextureOnSurfaceFilter::RequestData(vtkInformation* request,
46 vtkInformationVector** inputVector,
47 vtkInformationVector* outputVector)
48{
49 (void)request;
50
51 if (this->GetTexture() == nullptr)
52 {
53 ERR("VtkTextureOnSurfaceFilter::RequestData() - No texture specified.");
54 return 0;
55 }
56
57 vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
58 vtkPolyData* input =
59 vtkPolyData::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT()));
60
61 int dims[3];
62 this->GetTexture()->GetInput()->GetDimensions(dims);
63 std::size_t imgWidth = dims[0];
64 std::size_t imgHeight = dims[1];
65
66 std::pair<int, int> min(static_cast<int>(_origin.first),
67 static_cast<int>(_origin.second));
68 std::pair<int, int> max(
69 static_cast<int>(_origin.first + (imgWidth * _scalingFactor)),
70 static_cast<int>(_origin.second + (imgHeight * _scalingFactor)));
71
72 // calculate texture coordinates
73 vtkPoints* points = input->GetPoints();
74 vtkSmartPointer<vtkFloatArray> textureCoordinates =
75 vtkSmartPointer<vtkFloatArray>::New();
76 textureCoordinates->SetNumberOfComponents(2);
77 std::size_t nPoints = points->GetNumberOfPoints();
78 textureCoordinates->SetNumberOfTuples(nPoints);
79 textureCoordinates->SetName("textureCoords");
80 /* // adaptation for netcdf-curtain for TERENO Demo
81 double dist(0.0);
82 for (std::size_t i = 0; i < nPoints; i++)
83 {
84 double coords[3];
85 if ((i==0) || (i==173))
86 {
87 if (i==0) dist=0;
88 }
89 else
90 {
91 points->GetPoint(i-1, coords);
92 GeoLib::Point* pnt = new GeoLib::Point(coords);
93 points->GetPoint(i, coords);
94 GeoLib::Point* pnt2 = new GeoLib::Point(coords);
95 if (i<173)
96 dist += std::sqrt(MathLib::sqrDist(pnt, pnt2));
97 else
98 dist -= std::sqrt(MathLib::sqrDist(pnt, pnt2));
99 }
100 points->GetPoint(i, coords);
101 double x = MathLib::normalize(0, 8404, dist);
102 double z = MathLib::normalize(-79.5, 1.5, coords[2]);
103 float newcoords[2] = {x, z};
104 textureCoordinates->InsertNextTuple(newcoords);
105 }
106 */
107
108 int const range[2] = {max.first - min.first, max.second - min.second};
109 // Scale values relative to the range.
110
111 double coords[3];
112 for (std::size_t i = 0; i < nPoints; i++)
113 {
114 points->GetPoint(i, coords);
115
116 textureCoordinates->SetTuple2(i,
117 (coords[0] - min.first) / range[0],
118 (coords[1] - min.second) / range[1]);
119 }
120
121 // put it all together
122 vtkInformation* outInfo = outputVector->GetInformationObject(0);
123 vtkPolyData* output =
124 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
125 output->CopyStructure(input);
126 output->GetPointData()->PassData(input->GetPointData());
127 output->GetCellData()->PassData(input->GetCellData());
128 output->GetPointData()->SetTCoords(textureCoordinates);
129 output->Squeeze();
130
131 return 1;
132}
133
134void VtkTextureOnSurfaceFilter::SetRaster(vtkImageAlgorithm* img)
135{
136 double range[2];
137 img->Update();
138 img->GetOutput()->GetPointData()->GetScalars()->GetRange(range);
139 vtkSmartPointer<vtkImageShiftScale> scale =
140 vtkSmartPointer<vtkImageShiftScale>::New();
141 scale->SetInputConnection(img->GetOutputPort());
142 scale->SetShift(-range[0]);
143 scale->SetScale(255.0 / (range[1] - range[0]));
144 scale->SetOutputScalarTypeToUnsignedChar(); // Comment this out to get
145 // colored grayscale textures
146 scale->Update();
147
148 vtkTexture* texture = vtkTexture::New();
149 texture->InterpolateOff();
150 texture->RepeatOff();
151 // texture->EdgeClampOn(); // does not work
152 texture->SetInputData(scale->GetOutput());
153 this->SetTexture(texture);
154
155 double* origin = img->GetOutput()->GetOrigin();
156 _origin = {origin[0], origin[1]};
157 _scalingFactor = img->GetOutput()->GetSpacing()[0];
158}
159
160void VtkTextureOnSurfaceFilter::SetUserProperty(QString name, QVariant value)
161{
163}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
vtkStandardNewMacro(VtkTextureOnSurfaceFilter)
Definition of the VtkTextureOnSurfaceFilter class.
Definition of the VtkVisHelper class.
vtkTexture * GetTexture()
Returns a texture (if one has been assigned).
void SetTexture(vtkTexture *t)
Sets a texture for the VtkVisPipelineItem.
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.
Filter class for assigning a texture to a surface.
std::pair< double, double > _origin
~VtkTextureOnSurfaceFilter() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void PrintSelf(ostream &os, vtkIndent indent) override
Prints the object data to an output stream.
void SetUserProperty(QString name, QVariant value) override
Sets a user property. This should be implemented by subclasses.
void SetRaster(vtkImageAlgorithm *img)
Sets the raster/image to be used as a texture map.