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

Static Public Member Functions

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

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

◆ ~VtkTextureOnSurfaceFilter()

VtkTextureOnSurfaceFilter::~VtkTextureOnSurfaceFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static VtkTextureOnSurfaceFilter* VtkTextureOnSurfaceFilter::New ( )
static

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

Referenced by VtkCompositeTextureOnSurfaceFilter::init().

◆ PrintSelf()

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

Prints the object data to an output stream.

Definition at line 40 of file VtkTextureOnSurfaceFilter.cpp.

41 {
42  this->Superclass::PrintSelf(os, indent);
43 }

◆ 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 45 of file VtkTextureOnSurfaceFilter.cpp.

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 }
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
vtkTexture * GetTexture()
Returns a texture (if one has been assigned).
std::pair< double, double > _origin

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 134 of file VtkTextureOnSurfaceFilter.cpp.

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

References _origin, _scalingFactor, MathLib::LinAlg::scale(), 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 160 of file VtkTextureOnSurfaceFilter.cpp.

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

References MaterialPropertyLib::name, and VtkAlgorithmProperties::SetUserProperty().

◆ vtkTypeMacro()

VtkTextureOnSurfaceFilter::vtkTypeMacro ( VtkTextureOnSurfaceFilter  ,
vtkPolyDataAlgorithm   
)

Member Data Documentation

◆ _origin

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

Definition at line 64 of file VtkTextureOnSurfaceFilter.h.

Referenced by RequestData(), and SetRaster().

◆ _scalingFactor

double VtkTextureOnSurfaceFilter::_scalingFactor {0.}
private

Definition at line 65 of file VtkTextureOnSurfaceFilter.h.

Referenced by RequestData(), and SetRaster().


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