OGS
VtkTextureOnSurfaceFilter Class Reference

Detailed Description

Filter class for assigning a texture to a surface.

Use SetRaster() to define the texture that should be mapped on the object. The input of this class is a vtkPolyData object. It is important to call SetTexture() before calling SetInputConnection(). Texture coordinates will then be calculated automatically and the texture will also be saved in the VtkAlgorithmProperties object from which this class is derived (i.e. the texture can be returned by VtkAlgorithmProperties::GetTexture()).

For convenience this class also has a converter function ConvertImageToTexture() which uses a QImage as input.

Definition at line 26 of file VtkTextureOnSurfaceFilter.h.

#include <VtkTextureOnSurfaceFilter.h>

Inheritance diagram for VtkTextureOnSurfaceFilter:
[legend]
Collaboration diagram for VtkTextureOnSurfaceFilter:
[legend]

Public Member Functions

 vtkTypeMacro (VtkTextureOnSurfaceFilter, vtkPolyDataAlgorithm)
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints the object data to an output stream.
void SetRaster (vtkImageAlgorithm *img)
 Sets the raster/image to be used as a texture map.
void SetUserProperty (QString name, QVariant value) override
 Sets a user property. This should be implemented by subclasses.
Public Member Functions inherited from VtkAlgorithmProperties
 VtkAlgorithmProperties (QObject *parent=nullptr)
 Constructor (sets default values)
 ~VtkAlgorithmProperties () override
vtkProperty * GetProperties () const
 Returns the vtk properties.
vtkTexture * GetTexture ()
 Returns a texture (if one has been assigned).
void SetTexture (vtkTexture *t)
 Sets a texture for the VtkVisPipelineItem.
vtkLookupTable * GetLookupTable (const QString &array_name)
 Returns the colour lookup table (if one has been assigned).
void RemoveLookupTable (const QString &array_name)
 Removes the lookup table for the given scalar.
void SetLookUpTable (const QString &array_name, vtkLookupTable *lut)
 Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
void SetLookUpTable (const QString &array_name, const QString &filename)
 Loads a predefined color lookup table from a file for the specified scalar array.
bool GetScalarVisibility () const
 Returns the scalar visibility.
void SetScalarVisibility (bool on)
 Sets the scalar visibility.
QString GetName () const
 Returns the name. This is set to the file path if it is a source algorithm.
void SetName (QString name)
 Sets the name.
bool IsRemovable () const
 Is this algorithm removable from the pipeline (view).
QMap< QString, QVariant > * GetAlgorithmUserProperties () const
 Returns a map of user properties.
QMap< QString, QList< QVariant > > * GetAlgorithmUserVectorProperties () const
 Returns a map of vector user properties.
QVariant GetUserProperty (QString name) const
 Returns the value of a user property.
virtual void SetUserVectorProperty (QString name, QList< QVariant > values)
 Sets a vector user property. This should be implemented by subclasses.
QList< QVariant > GetUserVectorProperty (QString name) const
 Returns a list of values of a vector user property.
void SetActiveAttribute (QString name)
 Set the active attribute.
QString GetActiveAttribute () const
 Returns the desired active attribute.

Static Public Member Functions

static VtkTextureOnSurfaceFilterNew ()
 Create new objects with New() because of VTKs object reference counting.

Protected Member Functions

 VtkTextureOnSurfaceFilter ()
 ~VtkTextureOnSurfaceFilter () override
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override

Private Attributes

std::pair< double, double > _origin {0, 0}
double _scalingFactor {0.}

Additional Inherited Members

Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
Protected Attributes inherited from VtkAlgorithmProperties
vtkProperty * _property
vtkTexture * _texture
bool _scalarVisibility
std::map< QString, vtkLookupTable * > _lut
QString _name
QString _activeAttributeName
bool _removable
QMap< QString, QVariant > * _algorithmUserProperties
QMap< QString, QList< QVariant > > * _algorithmUserVectorProperties

Constructor & Destructor Documentation

◆ VtkTextureOnSurfaceFilter()

VtkTextureOnSurfaceFilter::VtkTextureOnSurfaceFilter ( )
protecteddefault

Referenced by New(), and vtkTypeMacro().

◆ ~VtkTextureOnSurfaceFilter()

VtkTextureOnSurfaceFilter::~VtkTextureOnSurfaceFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

VtkTextureOnSurfaceFilter * VtkTextureOnSurfaceFilter::New ( )
static

Create new objects with New() because of VTKs object reference counting.

References VtkTextureOnSurfaceFilter().

Referenced by VtkCompositeTextureOnSurfaceFilter::init().

◆ PrintSelf()

void VtkTextureOnSurfaceFilter::PrintSelf ( ostream & os,
vtkIndent indent )
override

Prints the object data to an output stream.

Definition at line 29 of file VtkTextureOnSurfaceFilter.cpp.

30{
31 this->Superclass::PrintSelf(os, indent);
32}

◆ RequestData()

int VtkTextureOnSurfaceFilter::RequestData ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Copies the input data and calculates texture coordinates (this requires that SetRaster() has been called before this method is executed.

Definition at line 34 of file VtkTextureOnSurfaceFilter.cpp.

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}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
vtkTexture * GetTexture()
Returns a texture (if one has been assigned).
std::pair< double, double > _origin
constexpr ranges::views::view_closure coords
Definition Mesh.h:223

References _origin, _scalingFactor, ERR(), and VtkAlgorithmProperties::GetTexture().

◆ SetRaster()

void VtkTextureOnSurfaceFilter::SetRaster ( vtkImageAlgorithm * img)

Sets the raster/image to be used as a texture map.

Definition at line 123 of file VtkTextureOnSurfaceFilter.cpp.

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}
void SetTexture(vtkTexture *t)
Sets a texture for the VtkVisPipelineItem.
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:37

References _origin, _scalingFactor, and VtkAlgorithmProperties::SetTexture().

Referenced by VtkCompositeTextureOnSurfaceFilter::init().

◆ SetUserProperty()

void VtkTextureOnSurfaceFilter::SetUserProperty ( QString name,
QVariant value )
overridevirtual

Sets a user property. This should be implemented by subclasses.

Reimplemented from VtkAlgorithmProperties.

Definition at line 149 of file VtkTextureOnSurfaceFilter.cpp.

150{
152}
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.

References VtkAlgorithmProperties::SetUserProperty().

◆ vtkTypeMacro()

VtkTextureOnSurfaceFilter::vtkTypeMacro ( VtkTextureOnSurfaceFilter ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ _origin

std::pair<double, double> VtkTextureOnSurfaceFilter::_origin {0, 0}
private

Definition at line 53 of file VtkTextureOnSurfaceFilter.h.

53{0, 0};

Referenced by RequestData(), and SetRaster().

◆ _scalingFactor

double VtkTextureOnSurfaceFilter::_scalingFactor {0.}
private

Definition at line 54 of file VtkTextureOnSurfaceFilter.h.

54{0.};

Referenced by RequestData(), and SetRaster().


The documentation for this class was generated from the following files: