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