OGS
VtkColorByHeightFilter Class Reference

Detailed Description

VTK filter object for colouring vtkPolyData objects based on z-coordinates.

This filter class is basically a container for a ColorLookupTable. In fact, you can get the underlying ColorLookupTable using the method GetColorLookupTable(). Using this method allows the user to set a number of properties on that lookup table such as interpolation method, the range of values over which the lookup table is calculated and so on. If no range boundaries are explicitly set, the minimum and maximum height value will be calculated from the data and set as minimum and maximum value for the lookup table. ColorLookupTable must be deleted manually.

See also
VtkCompositeColorByHeightFilter::init() for example usage.

Definition at line 37 of file VtkColorByHeightFilter.h.

#include <VtkColorByHeightFilter.h>

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

Public Member Functions

 vtkTypeMacro (VtkColorByHeightFilter, vtkPolyDataAlgorithm)
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints the mesh data to an output stream.
 
 vtkGetObjectMacro (ColorLookupTable, VtkColorLookupTable)
 Returns the underlying colour look up table object.
 
vtkMTimeType GetMTime () override
 This filter gets updated when the color look-up table was modified.
 
void SetUserProperty (QString name, QVariant value) override
 Sets user properties.
 
void SetTableRange (double min, double max)
 Sets the boundaries for the color look-up table.
 
void SetTableRangeScaling (double scale)
 Sets the scaling of the color look-up table boundaries. This is used in VtkVisTabWidget when a parent filter is scaled.
 
- 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 VtkColorByHeightFilterNew ()
 Create new objects with New() because of VTKs object reference counting.
 

Protected Member Functions

 VtkColorByHeightFilter ()
 
 ~VtkColorByHeightFilter () override
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 The filter logic.
 
VtkColorLookupTableBuildColorTable ()
 Calculates the color lookup table based on set parameters.
 

Protected Attributes

VtkColorLookupTableColorLookupTable
 
double _tableRange [2]
 
double _tableRangeScaling {1.0}
 
- 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
 

Additional Inherited Members

- Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
 

Constructor & Destructor Documentation

◆ VtkColorByHeightFilter()

VtkColorByHeightFilter::VtkColorByHeightFilter ( )
protected

Definition at line 34 of file VtkColorByHeightFilter.cpp.

36{
37}
VtkColorLookupTable * ColorLookupTable
static VtkColorLookupTable * New()
Create new objects with New() because of VTKs object reference counting.

◆ ~VtkColorByHeightFilter()

VtkColorByHeightFilter::~VtkColorByHeightFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ BuildColorTable()

VtkColorLookupTable * VtkColorByHeightFilter::BuildColorTable ( )
protected

Calculates the color lookup table based on set parameters.

◆ GetMTime()

vtkMTimeType VtkColorByHeightFilter::GetMTime ( )
override

This filter gets updated when the color look-up table was modified.

Definition at line 53 of file VtkColorByHeightFilter.cpp.

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}

References ColorLookupTable.

◆ New()

static VtkColorByHeightFilter * VtkColorByHeightFilter::New ( )
static

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

Referenced by VtkCompositeColorByHeightFilter::init().

◆ PrintSelf()

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

Prints the mesh data to an output stream.

Definition at line 41 of file VtkColorByHeightFilter.cpp.

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}
DataHolderLib::LUTType getInterpolationType() const
Returns the type of interpolation used.

References ColorLookupTable, and VtkColorLookupTable::getInterpolationType().

◆ RequestData()

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

The filter logic.

Definition at line 63 of file VtkColorByHeightFilter.cpp.

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}

◆ SetTableRange()

void VtkColorByHeightFilter::SetTableRange ( double min,
double max )

Sets the boundaries for the color look-up table.

Definition at line 98 of file VtkColorByHeightFilter.cpp.

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}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45

References _tableRange, ColorLookupTable, and ERR().

◆ SetTableRangeScaling()

void VtkColorByHeightFilter::SetTableRangeScaling ( double scale)

Sets the scaling of the color look-up table boundaries. This is used in VtkVisTabWidget when a parent filter is scaled.

Definition at line 113 of file VtkColorByHeightFilter.cpp.

114{
116 this->ColorLookupTable->SetTableRange(this->_tableRange[0] * scale,
117 this->_tableRange[1] * scale);
118}
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:44

References _tableRange, _tableRangeScaling, and ColorLookupTable.

◆ SetUserProperty()

void VtkColorByHeightFilter::SetUserProperty ( QString name,
QVariant value )
inlineoverridevirtual

Sets user properties.

Reimplemented from VtkAlgorithmProperties.

Definition at line 55 of file VtkColorByHeightFilter.h.

56 {
57 Q_UNUSED(name);
58 Q_UNUSED(value);
59 }

◆ vtkGetObjectMacro()

VtkColorByHeightFilter::vtkGetObjectMacro ( ColorLookupTable ,
VtkColorLookupTable  )

Returns the underlying colour look up table object.

◆ vtkTypeMacro()

VtkColorByHeightFilter::vtkTypeMacro ( VtkColorByHeightFilter ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ _tableRange

double VtkColorByHeightFilter::_tableRange[2]
protected

Definition at line 82 of file VtkColorByHeightFilter.h.

Referenced by SetTableRange(), and SetTableRangeScaling().

◆ _tableRangeScaling

double VtkColorByHeightFilter::_tableRangeScaling {1.0}
protected

Definition at line 83 of file VtkColorByHeightFilter.h.

83{1.0};

Referenced by SetTableRangeScaling().

◆ ColorLookupTable

VtkColorLookupTable* VtkColorByHeightFilter::ColorLookupTable
protected

Definition at line 80 of file VtkColorByHeightFilter.h.

Referenced by GetMTime(), PrintSelf(), SetTableRange(), and SetTableRangeScaling().


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