OGS
VtkPointsSource Class Reference

Detailed Description

VtkPointsSource is a VTK source object for the visualization of point data. As a vtkPolyDataAlgorithm it outputs polygonal data.

Definition at line 27 of file VtkPointsSource.h.

#include <VtkPointsSource.h>

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

Public Member Functions

 vtkTypeMacro (VtkPointsSource, vtkPolyDataAlgorithm)
 
void setPoints (const std::vector< GeoLib::Point * > *points)
 Sets the points as a vector.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints its data on a stream.
 
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 VtkPointsSourceNew ()
 Create new objects with New() because of VTKs object reference counting.
 

Protected Member Functions

 VtkPointsSource ()
 
 ~VtkPointsSource () override=default
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 Computes the polygonal data object.
 
int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 

Protected Attributes

const std::vector< GeoLib::Point * > * _points {nullptr}
 The points to visualize.
 
- 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

◆ VtkPointsSource()

VtkPointsSource::VtkPointsSource ( )
protected

Definition at line 34 of file VtkPointsSource.cpp.

35{
36 _removable = false; // From VtkAlgorithmProperties
37 this->SetNumberOfInputPorts(0);
38
40 GetProperties()->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
41}
vtkProperty * GetProperties() const
Returns the vtk properties.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24

References VtkAlgorithmProperties::_removable, VtkAlgorithmProperties::GetProperties(), and DataHolderLib::getRandomColor().

◆ ~VtkPointsSource()

VtkPointsSource::~VtkPointsSource ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static VtkPointsSource * VtkPointsSource::New ( )
static

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

◆ PrintSelf()

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

Prints its data on a stream.

Definition at line 43 of file VtkPointsSource.cpp.

44{
45 this->Superclass::PrintSelf(os, indent);
46
47 if (_points->empty())
48 {
49 return;
50 }
51
52 os << indent << "== VtkPointsSource =="
53 << "\n";
54
55 int i = 0;
56 for (auto point : *_points)
57 {
58 os << indent << "Point " << i << " (" << (*point)[0] << ", "
59 << (*point)[1] << ", " << (*point)[2] << ")\n";
60 i++;
61 }
62}
const std::vector< GeoLib::Point * > * _points
The points to visualize.

References _points.

◆ RequestData()

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

Computes the polygonal data object.

Definition at line 64 of file VtkPointsSource.cpp.

67{
68 (void)request;
69 (void)inputVector;
70
71 if (!_points)
72 {
73 return 0;
74 }
75 int numPoints = _points->size();
76 if (numPoints == 0)
77 {
78 ERR("VtkPointsSource::RequestData(): Size of point vector is 0");
79 return 0;
80 }
81
82 vtkSmartPointer<vtkInformation> outInfo =
83 outputVector->GetInformationObject(0);
84 vtkSmartPointer<vtkPolyData> output =
85 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
86
87 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
88 vtkSmartPointer<vtkCellArray> newVerts =
89 vtkSmartPointer<vtkCellArray>::New();
90 newPoints->SetNumberOfPoints(numPoints);
91 newVerts->Allocate(numPoints);
92
93 vtkSmartPointer<vtkIntArray> pointIDs = vtkSmartPointer<vtkIntArray>::New();
94 pointIDs->SetNumberOfComponents(1);
95 pointIDs->SetNumberOfValues(numPoints);
96 pointIDs->SetName("PointIDs");
97
98 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
99 0)
100 {
101 return 1;
102 }
103
104 // Generate points and vertices
105 unsigned i = 0;
106 for (auto point : *_points)
107 {
108 double coords[3] = {(*point)[0], (*point)[1], (*point)[2]};
109 newPoints->SetPoint(i, coords);
110 newVerts->InsertNextCell(1);
111 newVerts->InsertCellPoint(i);
112
113 pointIDs->SetValue(i, i);
114 i++;
115 }
116
117 output->SetPoints(newPoints);
118 output->SetVerts(newVerts);
119 output->GetCellData()->AddArray(pointIDs);
120 output->GetCellData()->SetActiveAttribute("PointIDs",
121 vtkDataSetAttributes::SCALARS);
122
123 return 1;
124}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
constexpr ranges::views::view_closure coords
Definition Mesh.h:232

References _points, and ERR().

◆ RequestInformation()

int VtkPointsSource::RequestInformation ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Definition at line 126 of file VtkPointsSource.cpp.

129{
130 return 1;
131}

◆ setPoints()

void VtkPointsSource::setPoints ( const std::vector< GeoLib::Point * > * points)
inline

Sets the points as a vector.

Definition at line 36 of file VtkPointsSource.h.

36{ _points = points; }

References _points.

◆ SetUserProperty()

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

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

Reimplemented from VtkAlgorithmProperties.

Definition at line 133 of file VtkPointsSource.cpp.

134{
135 Q_UNUSED(name);
136 Q_UNUSED(value);
137}

◆ vtkTypeMacro()

VtkPointsSource::vtkTypeMacro ( VtkPointsSource ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ _points

const std::vector<GeoLib::Point*>* VtkPointsSource::_points {nullptr}
protected

The points to visualize.

Definition at line 57 of file VtkPointsSource.h.

57{nullptr};

Referenced by PrintSelf(), RequestData(), and setPoints().


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